libMesh::LaspackVector< T > Class Template Reference

#include <laspack_matrix.h>

Inheritance diagram for libMesh::LaspackVector< T >:

Public Member Functions

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

Static Public Member Functions

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

Protected Types

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

Protected Member Functions

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

Protected Attributes

bool _is_closed
 
bool _is_initialized
 
ParallelType _type
 
const Parallel::Communicator_communicator
 

Static Protected Attributes

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

Private Attributes

QVector _vec
 

Friends

class LaspackLinearSolver< T >
 

Detailed Description

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

This class provides a nice interface to the Laspack C-based data structures for serial vectors. All overridden virtual functions are documented in numeric_vector.h.

Author
Benjamin S. Kirk
Date
2002

Definition at line 42 of file laspack_matrix.h.

Member Typedef Documentation

◆ Counts

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

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

Definition at line 117 of file reference_counter.h.

Constructor & Destructor Documentation

◆ LaspackVector() [1/6]

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

Dummy-Constructor. Dimension=0

Definition at line 248 of file laspack_vector.h.

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

250  : NumericVector<T>(comm, ptype)
251 {
252  this->_type = ptype;
253 }
const Parallel::Communicator & comm() const

◆ LaspackVector() [2/6]

template<typename T >
libMesh::LaspackVector< T >::LaspackVector ( 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 259 of file laspack_vector.h.

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

262  : NumericVector<T>(comm, ptype)
263 {
264  this->init(n, n, false, ptype);
265 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
const Parallel::Communicator & comm() const

◆ LaspackVector() [3/6]

template<typename T >
libMesh::LaspackVector< T >::LaspackVector ( 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 271 of file laspack_vector.h.

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

275  : NumericVector<T>(comm, ptype)
276 {
277  this->init(n, n_local, false, ptype);
278 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
const Parallel::Communicator & comm() const

◆ LaspackVector() [4/6]

template<typename T >
libMesh::LaspackVector< T >::LaspackVector ( 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 284 of file laspack_vector.h.

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

289  : NumericVector<T>(comm, ptype)
290 {
291  this->init(N, n_local, ghost, false, ptype);
292 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
const Parallel::Communicator & comm() const

◆ LaspackVector() [5/6]

template<typename T>
libMesh::LaspackVector< T >::LaspackVector ( LaspackVector< T > &&  )
delete

This class manages a C-style struct (QVector) manually, so we don't want to allow any automatic copy/move functions to be generated, and we can't default the destructor.

◆ LaspackVector() [6/6]

template<typename T>
libMesh::LaspackVector< T >::LaspackVector ( const LaspackVector< T > &  )
delete

◆ ~LaspackVector()

template<typename T >
libMesh::LaspackVector< T >::~LaspackVector ( )
inlinevirtual

Definition at line 298 of file laspack_vector.h.

299 {
300  this->clear ();
301 }
virtual void clear() override

Member Function Documentation

◆ abs()

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

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

Implements libMesh::NumericVector< T >.

Definition at line 237 of file laspack_vector.C.

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

238 {
239  libmesh_assert (this->initialized());
240 
241  const numeric_index_type n = this->size();
242 
243  for (numeric_index_type i=0; i!=n; ++i)
244  this->set(i,std::abs((*this)(i)));
245 }
double abs(double a)
virtual numeric_index_type size() const override
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ add() [1/4]

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

Adds value to each entry of the vector.

Implements libMesh::NumericVector< T >.

Definition at line 510 of file laspack_vector.h.

References libMesh::initialized(), and value.

511 {
512  libmesh_assert (this->initialized());
513  libmesh_assert_less (i, this->size());
514 
515  V_AddCmp (&_vec, i+1, value);
516 
517 #ifndef NDEBUG
518  this->_is_closed = false;
519 #endif
520 }
virtual numeric_index_type size() const override
virtual bool initialized() const
static const bool value
Definition: xdr_io.C:109

◆ add() [2/4]

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

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

Implements libMesh::NumericVector< T >.

Definition at line 156 of file laspack_vector.C.

157 {
158  const numeric_index_type n = this->size();
159 
160  for (numeric_index_type i=0; i<n; i++)
161  this->add (i, v);
162 
163 #ifndef NDEBUG
164  this->_is_closed = false;
165 #endif
166 }
virtual numeric_index_type size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual void add(const numeric_index_type i, const T value) override

◆ add() [3/4]

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

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

Implements libMesh::NumericVector< T >.

Definition at line 172 of file laspack_vector.C.

173 {
174  this->add (1., v);
175 }
virtual void add(const numeric_index_type i, const T value) override

◆ add() [4/4]

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

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

Implements libMesh::NumericVector< T >.

Definition at line 180 of file laspack_vector.C.

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

181 {
182  // Make sure the vector passed in is really a LaspackVector
183  const LaspackVector * v = cast_ptr<const LaspackVector *>(&v_in);
184 
185 #ifndef NDEBUG
186  const bool was_closed = this->_is_closed;
187 #endif
188 
189  libmesh_assert(v);
190  libmesh_assert_equal_to (this->size(), v->size());
191 
192  for (numeric_index_type i=0; i<v->size(); i++)
193  this->add (i, a*(*v)(i));
194 
195 #ifndef NDEBUG
196  this->_is_closed = was_closed;
197 #endif
198 }
virtual numeric_index_type size() const override
LaspackVector(const Parallel::Communicator &comm, const ParallelType=AUTOMATIC)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual void add(const numeric_index_type i, const T value) override

◆ add_vector() [1/6]

template<typename T >
void libMesh::LaspackVector< T >::add_vector ( const NumericVector< T > &  v,
const SparseMatrix< T > &  A 
)
overridevirtual

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

Implements libMesh::NumericVector< T >.

Definition at line 203 of file laspack_vector.C.

References libMesh::LaspackVector< T >::_vec.

205 {
206  // Make sure the data passed in are really in Laspack types
207  const LaspackVector<T> * vec = cast_ptr<const LaspackVector<T> *>(&vec_in);
208  const LaspackMatrix<T> * mat = cast_ptr<const LaspackMatrix<T> *>(&mat_in);
209 
210  libmesh_assert(vec);
211  libmesh_assert(mat);
212 
213  // += mat*vec
214  AddAsgn_VV (&_vec, Mul_QV(const_cast<QMatrix*>(&mat->_QMat),
215  const_cast<QVector*>(&vec->_vec)));
216 }

◆ add_vector() [2/6]

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

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

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

Definition at line 362 of file numeric_vector.C.

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

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

◆ add_vector() [3/6]

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

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

Definition at line 835 of file numeric_vector.h.

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

◆ add_vector() [4/6]

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

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

Definition at line 372 of file numeric_vector.C.

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

◆ add_vector() [5/6]

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

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

Definition at line 847 of file numeric_vector.h.

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

◆ add_vector() [6/6]

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

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

Definition at line 384 of file numeric_vector.C.

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

◆ add_vector_transpose()

template<typename T >
void libMesh::LaspackVector< T >::add_vector_transpose ( const NumericVector< T > &  v,
const SparseMatrix< T > &  A 
)
overridevirtual

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

Implements libMesh::NumericVector< T >.

Definition at line 220 of file laspack_vector.C.

222 {
223  libmesh_not_implemented();
224 }

◆ build()

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

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

Definition at line 48 of file numeric_vector.C.

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

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

◆ clear()

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

Restores the NumericVector<T> to a pristine state.

Reimplemented from libMesh::NumericVector< T >.

Definition at line 394 of file laspack_vector.h.

References libMesh::libMeshPrivateData::_is_initialized, and libMesh::initialized().

395 {
396  if (this->initialized())
397  {
398  V_Destr (&_vec);
399  }
400 
401  this->_is_initialized = false;
402 #ifndef NDEBUG
403  this->_is_closed = false;
404 #endif
405 }
virtual bool initialized() const

◆ clone()

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

Implements libMesh::NumericVector< T >.

Definition at line 435 of file laspack_vector.h.

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

436 {
437  NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
438 
439  cloned_vector->init(*this, true);
440 
441  *cloned_vector = *this;
442 
443  return std::unique_ptr<NumericVector<T>>(cloned_vector);
444 }
const Parallel::Communicator & comm() const

◆ close()

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

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

Implements libMesh::NumericVector< T >.

Definition at line 381 of file laspack_vector.h.

References libMesh::initialized().

382 {
383  libmesh_assert (this->initialized());
384 
385 #ifndef NDEBUG
386  this->_is_closed = true;
387 #endif
388 }
virtual bool initialized() const

◆ closed()

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

◆ comm()

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

Definition at line 89 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

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

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

◆ compare()

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

Definition at line 104 of file numeric_vector.C.

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

◆ conjugate()

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

Negates the imaginary component of each entry in the vector.

Implements libMesh::NumericVector< T >.

Definition at line 142 of file laspack_vector.C.

References libMesh::libmesh_conj().

143 {
144  const numeric_index_type n = this->size();
145 
146  for (numeric_index_type i=0; i<n; i++)
147  {
148  T v = (*this)(i);
149 
150  this->set(i, libmesh_conj(v) );
151  }
152 }
T libmesh_conj(T a)
virtual numeric_index_type size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ create_subvector()

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

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

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

Definition at line 692 of file numeric_vector.h.

694  {
695  libmesh_not_implemented();
696  }

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

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

◆ dot()

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

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

Implements libMesh::NumericVector< T >.

Definition at line 248 of file laspack_vector.C.

References libMesh::LaspackVector< T >::_vec, and libMesh::initialized().

249 {
250  libmesh_assert (this->initialized());
251 
252  // Make sure the NumericVector passed in is really a LaspackVector
253  const LaspackVector<T> * v = cast_ptr<const LaspackVector<T> *>(&v_in);
254  libmesh_assert(v);
255 
256  return Mul_VV (const_cast<QVector*>(&(this->_vec)),
257  const_cast<QVector*>(&(v->_vec)));
258 }
virtual bool initialized() const

◆ el()

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

Definition at line 346 of file numeric_vector.h.

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

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

◆ first_local_index()

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

Implements libMesh::NumericVector< T >.

Definition at line 472 of file laspack_vector.h.

References libMesh::initialized().

473 {
474  libmesh_assert (this->initialized());
475 
476  return 0;
477 }
virtual bool initialized() const

◆ get() [1/2]

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

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

Reimplemented in libMesh::PetscVector< T >.

Definition at line 806 of file numeric_vector.h.

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

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

◆ get() [2/2]

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

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

Definition at line 820 of file numeric_vector.h.

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

◆ get_info()

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

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

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

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

◆ global_relative_compare()

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

Definition at line 169 of file numeric_vector.C.

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

◆ increment_constructor_count()

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

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

Definition at line 181 of file reference_counter.h.

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

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

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

◆ increment_destructor_count()

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

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

Definition at line 194 of file reference_counter.h.

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

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

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

◆ init() [1/4]

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

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

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

Implements libMesh::NumericVector< T >.

Referenced by libMesh::LaspackVector< T >::LaspackVector().

◆ init() [2/4]

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

Call init() with n_local = N.

Implements libMesh::NumericVector< T >.

Definition at line 346 of file laspack_vector.h.

References libMesh::TriangleWrapper::init().

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

◆ init() [3/4]

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

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

Implements libMesh::NumericVector< T >.

◆ init() [4/4]

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

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

Implements libMesh::NumericVector< T >.

Definition at line 371 of file laspack_vector.h.

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

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

◆ initialized()

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

◆ insert() [1/5]

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

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

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

Definition at line 82 of file numeric_vector.C.

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

◆ insert() [2/5]

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

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

Definition at line 859 of file numeric_vector.h.

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

◆ insert() [3/5]

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

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

Definition at line 92 of file numeric_vector.C.

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

◆ insert() [4/5]

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

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

Definition at line 871 of file numeric_vector.h.

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

◆ insert() [5/5]

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

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

Definition at line 883 of file numeric_vector.h.

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

◆ l1_norm()

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

Implements libMesh::NumericVector< T >.

Definition at line 54 of file laspack_vector.C.

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

55 {
56  libmesh_assert (this->closed());
57 
58  return static_cast<Real>(l1Norm_V(const_cast<QVector*>(&_vec)));
59 }
virtual bool closed() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ l2_norm()

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

Implements libMesh::NumericVector< T >.

Definition at line 64 of file laspack_vector.C.

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

65 {
66  libmesh_assert (this->closed());
67 
68  return static_cast<Real>(l2Norm_V(const_cast<QVector*>(&_vec)));
69 }
virtual bool closed() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ last_local_index()

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

Implements libMesh::NumericVector< T >.

Definition at line 483 of file laspack_vector.h.

References libMesh::initialized().

484 {
485  libmesh_assert (this->initialized());
486 
487  return this->size();
488 }
virtual numeric_index_type size() const override
virtual bool initialized() const

◆ linfty_norm()

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

Implements libMesh::NumericVector< T >.

Definition at line 74 of file laspack_vector.C.

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

75 {
76  libmesh_assert (this->closed());
77 
78  return static_cast<Real>(MaxNorm_V(const_cast<QVector*>(&_vec)));
79 }
virtual bool closed() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ local_relative_compare()

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

Definition at line 136 of file numeric_vector.C.

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

◆ local_size()

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

Implements libMesh::NumericVector< T >.

Definition at line 461 of file laspack_vector.h.

References libMesh::initialized().

462 {
463  libmesh_assert (this->initialized());
464 
465  return this->size();
466 }
virtual numeric_index_type size() const override
virtual bool initialized() const

◆ localize() [1/5]

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

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

Implements libMesh::NumericVector< T >.

Definition at line 395 of file laspack_vector.C.

397 {
398  v_local.resize(this->size());
399 
400  for (numeric_index_type i=0; i<v_local.size(); i++)
401  v_local[i] = (*this)(i);
402 }
virtual numeric_index_type size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ localize() [2/5]

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

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

Implements libMesh::NumericVector< T >.

Definition at line 335 of file laspack_vector.C.

336 {
337  // Make sure the NumericVector passed in is really a LaspackVector
338  LaspackVector<T> * v_local =
339  cast_ptr<LaspackVector<T> *>(&v_local_in);
340 
341  libmesh_assert(v_local);
342 
343  *v_local = *this;
344 }

◆ localize() [3/5]

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

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

Implements libMesh::NumericVector< T >.

◆ localize() [4/5]

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

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

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

Example:

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

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

Implements libMesh::NumericVector< T >.

Definition at line 365 of file laspack_vector.C.

367 {
368  // LaspackVectors are serial, so we can just copy values
369  v_local.resize(indices.size());
370 
371  for (numeric_index_type i=0; i<v_local.size(); i++)
372  v_local[i] = (*this)(indices[i]);
373 }
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ localize() [5/5]

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

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

Implements libMesh::NumericVector< T >.

◆ localize_to_one()

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

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

Implements libMesh::NumericVector< T >.

Definition at line 407 of file laspack_vector.C.

409 {
410  libmesh_assert_equal_to (pid, 0);
411 
412  this->localize (v_local);
413 }
virtual void localize(std::vector< T > &v_local) const override

◆ max()

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

Implements libMesh::NumericVector< T >.

Definition at line 427 of file laspack_vector.C.

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

428 {
429  libmesh_assert (this->initialized());
430  if (!this->size())
432 
433  Real the_max = libmesh_real((*this)(0));
434 
435  const numeric_index_type n = this->size();
436 
437  for (numeric_index_type i=1; i<n; i++)
438  the_max = std::max (the_max, libmesh_real((*this)(i)));
439 
440  return the_max;
441 }
T libmesh_real(T a)
virtual numeric_index_type size() const override
virtual bool initialized() const
long double max(long double a, double b)
dof_id_type numeric_index_type
Definition: id_types.h:92
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ min()

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

Implements libMesh::NumericVector< T >.

Definition at line 446 of file laspack_vector.C.

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

447 {
448  libmesh_assert (this->initialized());
449  if (!this->size())
451 
452  Real the_min = libmesh_real((*this)(0));
453 
454  const numeric_index_type n = this->size();
455 
456  for (numeric_index_type i=1; i<n; i++)
457  the_min = std::min (the_min, libmesh_real((*this)(i)));
458 
459  return the_min;
460 }
T libmesh_real(T a)
virtual numeric_index_type size() const override
virtual bool initialized() const
long double max(long double a, double b)
dof_id_type numeric_index_type
Definition: id_types.h:92
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
long double min(long double a, double b)

◆ n_objects()

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

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

Definition at line 83 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

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

◆ n_processors()

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

Definition at line 95 of file parallel_object.h.

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

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

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

◆ operator()()

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

Implements libMesh::NumericVector< T >.

Definition at line 526 of file laspack_vector.h.

References libMesh::initialized().

527 {
528  libmesh_assert (this->initialized());
529  libmesh_assert ( ((i >= this->first_local_index()) &&
530  (i < this->last_local_index())) );
531 
532 
533  return static_cast<T>(V_GetCmp(const_cast<QVector*>(&_vec), i+1));
534 }
virtual bool initialized() const
virtual numeric_index_type last_local_index() const override
virtual numeric_index_type first_local_index() const override

◆ operator*=()

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

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

Returns
A reference to *this.

Definition at line 391 of file numeric_vector.h.

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

◆ operator+=()

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

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

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 84 of file laspack_vector.C.

References libMesh::closed().

85 {
86  libmesh_assert (this->closed());
87 
88  this->add(1., v);
89 
90  return *this;
91 }
virtual void add(const numeric_index_type i, const T value) override
virtual bool closed() const

◆ operator-=()

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

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

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 97 of file laspack_vector.C.

References libMesh::closed().

98 {
99  libmesh_assert (this->closed());
100 
101  this->add(-1., v);
102 
103  return *this;
104 }
virtual void add(const numeric_index_type i, const T value) override
virtual bool closed() const

◆ operator/=() [1/2]

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

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

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 109 of file laspack_vector.C.

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

110 {
111  libmesh_assert_equal_to(size(), v.size());
112 
113  const numeric_index_type n = this->size();
114 
115  for (numeric_index_type i=0; i<n; i++)
116  this->set(i, (*this)(i) / v(i));
117 
118  return *this;
119 }
virtual numeric_index_type size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ operator/=() [2/2]

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

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

Returns
A reference to *this.

Definition at line 400 of file numeric_vector.h.

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

◆ operator=() [1/5]

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

Copy assignment operator. Calls Asgn_VV() to assign the contents of one vector to another.

Returns
A reference to *this as the derived type.

Definition at line 295 of file laspack_vector.C.

References libMesh::LaspackVector< T >::_vec, libMesh::NumericVector< T >::closed(), libMesh::initialized(), and libMesh::LaspackVector< T >::size().

296 {
297  libmesh_assert (this->initialized());
298  libmesh_assert (v.closed());
299  libmesh_assert_equal_to (this->size(), v.size());
300 
301  if (v.size() != 0)
302  Asgn_VV (const_cast<QVector*>(&_vec),
303  const_cast<QVector*>(&v._vec)
304  );
305 
306 #ifndef NDEBUG
307  this->_is_closed = true;
308 #endif
309 
310  return *this;
311 }
virtual numeric_index_type size() const override
virtual bool initialized() const

◆ operator=() [2/5]

template<typename T>
LaspackVector& libMesh::LaspackVector< T >::operator= ( LaspackVector< T > &&  )
delete

◆ operator=() [3/5]

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

Sets all entries of the vector to the value s.

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 264 of file laspack_vector.C.

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

265 {
266  libmesh_assert (this->initialized());
267  libmesh_assert (this->closed());
268 
269  V_SetAllCmp (&_vec, s);
270 
271  return *this;
272 }
virtual bool initialized() const
virtual bool closed() const

◆ operator=() [4/5]

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

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

Returns
A reference to *this as the base type.

Implements libMesh::NumericVector< T >.

Definition at line 278 of file laspack_vector.C.

279 {
280  // Make sure the NumericVector passed in is really a LaspackVector
281  const LaspackVector<T> * v =
282  cast_ptr<const LaspackVector<T> *>(&v_in);
283 
284  libmesh_assert(v);
285 
286  *this = *v;
287 
288  return *this;
289 }

◆ operator=() [5/5]

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

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

Returns
A reference to *this as the base type.

Case 1: The vector is the same size of The global vector. Only add the local components.

Implements libMesh::NumericVector< T >.

Definition at line 317 of file laspack_vector.C.

318 {
323  if (this->size() == v.size())
324  for (numeric_index_type i=0; i<v.size(); i++)
325  this->set (i, v[i]);
326 
327  else
328  libmesh_error_msg("this->size() = " << this->size() << " must be equal to v.size() = " << v.size());
329 
330  return *this;
331 }
virtual numeric_index_type size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ pointwise_mult()

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

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

Implements libMesh::NumericVector< T >.

Definition at line 418 of file laspack_vector.C.

420 {
421  libmesh_not_implemented();
422 }

◆ print() [1/2]

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

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

Definition at line 916 of file numeric_vector.h.

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

◆ print() [2/2]

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

Definition at line 898 of file numeric_vector.h.

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

◆ print_global() [1/2]

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

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

Definition at line 953 of file numeric_vector.h.

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

◆ print_global() [2/2]

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

Definition at line 931 of file numeric_vector.h.

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

◆ print_info()

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

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

Definition at line 87 of file reference_counter.C.

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

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

◆ print_matlab()

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

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

Reimplemented in libMesh::PetscVector< T >.

Definition at line 681 of file numeric_vector.h.

682  {
683  libmesh_not_implemented();
684  }

◆ processor_id()

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

Definition at line 101 of file parallel_object.h.

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

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

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

◆ reciprocal()

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

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

Implements libMesh::NumericVector< T >.

Definition at line 124 of file laspack_vector.C.

125 {
126  const numeric_index_type n = this->size();
127 
128  for (numeric_index_type i=0; i<n; i++)
129  {
130  T v = (*this)(i);
131 
132  // Don't divide by zero!
133  libmesh_assert_not_equal_to (v, T(0));
134 
135  this->set(i, 1. / v);
136  }
137 }
virtual numeric_index_type size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ scale()

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

Scale each element of the vector by the given factor.

Implements libMesh::NumericVector< T >.

Definition at line 229 of file laspack_vector.C.

References libMesh::initialized().

230 {
231  libmesh_assert (this->initialized());
232 
233  Asgn_VV(&_vec, Mul_SV (factor, &_vec));
234 }
virtual bool initialized() const

◆ set()

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

Sets v(i) = value.

Implements libMesh::NumericVector< T >.

Definition at line 494 of file laspack_vector.h.

References libMesh::initialized(), and value.

495 {
496  libmesh_assert (this->initialized());
497  libmesh_assert_less (i, this->size());
498 
499  V_SetCmp (&_vec, i+1, value);
500 
501 #ifndef NDEBUG
502  this->_is_closed = false;
503 #endif
504 }
virtual numeric_index_type size() const override
virtual bool initialized() const
static const bool value
Definition: xdr_io.C:109

◆ size()

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

Implements libMesh::NumericVector< T >.

Definition at line 450 of file laspack_vector.h.

References libMesh::initialized().

Referenced by libMesh::LaspackVector< T >::add(), and libMesh::LaspackVector< T >::operator=().

451 {
452  libmesh_assert (this->initialized());
453 
454  return static_cast<numeric_index_type>(V_GetDim(const_cast<QVector*>(&_vec)));
455 }
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ subset_l1_norm()

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

Definition at line 311 of file numeric_vector.C.

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

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

◆ subset_l2_norm()

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

Definition at line 326 of file numeric_vector.C.

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

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

◆ subset_linfty_norm()

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

Definition at line 341 of file numeric_vector.C.

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

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

◆ sum()

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

Implements libMesh::NumericVector< T >.

Definition at line 37 of file laspack_vector.C.

References libMesh::closed().

38 {
39  libmesh_assert (this->closed());
40 
41  T _sum = 0;
42 
43  const numeric_index_type n = this->size();
44 
45  for (numeric_index_type i=0; i!=n; ++i)
46  _sum += (*this)(i);
47 
48  return _sum;
49 }
virtual numeric_index_type size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual bool closed() const

◆ swap()

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

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

Reimplemented from libMesh::NumericVector< T >.

Definition at line 540 of file laspack_vector.h.

References swap().

541 {
542  LaspackVector<T> & v = cast_ref<LaspackVector<T> &>(other);
543 
544  // This is all grossly dependent on Laspack version...
545 
546  std::swap(_vec.Name, v._vec.Name);
547  std::swap(_vec.Dim, v._vec.Dim);
548  std::swap(_vec.Instance, v._vec.Instance);
549  std::swap(_vec.LockLevel, v._vec.LockLevel);
550  std::swap(_vec.Multipl, v._vec.Multipl);
551  std::swap(_vec.OwnData, v._vec.OwnData);
552 
553  // This should still be O(1), since _vec.Cmp is just a pointer to
554  // data on the heap
555 
556  std::swap(_vec.Cmp, v._vec.Cmp);
557 }
void swap(Iterator &lhs, Iterator &rhs)

◆ type() [1/2]

◆ type() [2/2]

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

Definition at line 159 of file numeric_vector.h.

159 { return _type; }

◆ zero()

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

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

Implements libMesh::NumericVector< T >.

Definition at line 410 of file laspack_vector.h.

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

Referenced by libMesh::LaspackLinearSolver< T >::adjoint_solve(), and libMesh::LaspackLinearSolver< T >::solve().

411 {
412  libmesh_assert (this->initialized());
413  libmesh_assert (this->closed());
414 
415  V_SetAllCmp (&_vec, 0.);
416 }
virtual bool initialized() const
virtual bool closed() const

◆ zero_clone()

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

Implements libMesh::NumericVector< T >.

Definition at line 422 of file laspack_vector.h.

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

423 {
424  NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
425 
426  cloned_vector->init(*this);
427 
428  return std::unique_ptr<NumericVector<T>>(cloned_vector);
429 }
const Parallel::Communicator & comm() const

Friends And Related Function Documentation

◆ LaspackLinearSolver< T >

template<typename T>
friend class LaspackLinearSolver< T >
friend

Make other Laspack datatypes friends

Definition at line 239 of file laspack_vector.h.

Member Data Documentation

◆ _communicator

◆ _counts

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

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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

Definition at line 141 of file reference_counter.h.

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

◆ _is_closed

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

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

Definition at line 712 of file numeric_vector.h.

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

◆ _is_initialized

◆ _mutex

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

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

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

Definition at line 130 of file reference_counter.h.

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

◆ _type

◆ _vec

template<typename T>
QVector libMesh::LaspackVector< T >::_vec
private

Actual Laspack vector datatype to hold vector entries

Definition at line 234 of file laspack_vector.h.

Referenced by libMesh::LaspackVector< T >::add_vector(), libMesh::LaspackVector< T >::dot(), and libMesh::LaspackVector< T >::operator=().


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