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

Static Public Member Functions

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

Protected Types

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

Protected Member Functions

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

Protected Attributes

bool _is_closed
 
bool _is_initialized
 
ParallelType _type
 
const Parallel::Communicator_communicator
 

Static Protected Attributes

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

Private Attributes

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

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

Constructor & Destructor Documentation

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

Dummy-Constructor. Dimension=0

Definition at line 244 of file laspack_vector.h.

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

246  : NumericVector<T>(comm, ptype)
247 {
248  this->_type = ptype;
249 }
const Parallel::Communicator & comm() const
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 255 of file laspack_vector.h.

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

258  : NumericVector<T>(comm, ptype)
259 {
260  this->init(n, n, false, ptype);
261 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
const Parallel::Communicator & comm() const
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 267 of file laspack_vector.h.

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

271  : NumericVector<T>(comm, ptype)
272 {
273  this->init(n, n_local, false, ptype);
274 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
const Parallel::Communicator & comm() const
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 280 of file laspack_vector.h.

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

285  : NumericVector<T>(comm, ptype)
286 {
287  this->init(N, n_local, ghost, false, ptype);
288 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
const Parallel::Communicator & comm() const

Member Function Documentation

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

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

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 bool initialized() const
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::add ( const numeric_index_type  i,
const T  value 
)
inlinevirtual

Adds value to each entry of the vector.

Implements libMesh::NumericVector< T >.

Definition at line 506 of file laspack_vector.h.

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

507 {
508  libmesh_assert (this->initialized());
509  libmesh_assert_less (i, this->size());
510 
511  V_AddCmp (&_vec, i+1, value);
512 
513 #ifndef NDEBUG
514  this->_is_closed = false;
515 #endif
516 }
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::add ( const T  s)
virtual

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 void add(const numeric_index_type i, const T value) libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::add ( const NumericVector< T > &  v)
virtual

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) libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::add ( const T  a,
const NumericVector< T > &  v 
)
virtual

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::libmesh_assert(), and 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 }
LaspackVector(const Parallel::Communicator &comm, const ParallelType=AUTOMATIC)
libmesh_assert(j)
virtual void add(const numeric_index_type i, const T value) libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::add_vector ( const NumericVector< T > &  v,
const SparseMatrix< T > &  A 
)
virtual

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, and libMesh::libmesh_assert().

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 }
libmesh_assert(j)
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 383 of file numeric_vector.C.

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

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

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 837 of file numeric_vector.h.

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

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

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 394 of file numeric_vector.C.

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

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

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 849 of file numeric_vector.h.

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

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

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 406 of file numeric_vector.C.

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

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

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 }
template<typename T >
UniquePtr< NumericVector< T > > libMesh::NumericVector< T >::build ( const Parallel::Communicator comm,
const SolverPackage  solver_package = libMesh::default_solver_package() 
)
staticinherited

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

Definition at line 46 of file numeric_vector.C.

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

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

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

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

Deprecated:
LIBMESH_DISABLE_COMMWORLD is now the default, use the build() method that takes a Parallel::Communicator instead.

Definition at line 84 of file numeric_vector.C.

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

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

Restores the NumericVector<T> to a pristine state.

Reimplemented from libMesh::NumericVector< T >.

Definition at line 390 of file laspack_vector.h.

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

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

391 {
392  if (this->initialized())
393  {
394  V_Destr (&_vec);
395  }
396 
397  this->_is_initialized = false;
398 #ifndef NDEBUG
399  this->_is_closed = false;
400 #endif
401 }
virtual bool initialized() const
template<typename T >
UniquePtr< NumericVector< T > > libMesh::LaspackVector< T >::clone ( ) const
inlinevirtual
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 431 of file laspack_vector.h.

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

432 {
433  NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
434 
435  cloned_vector->init(*this, true);
436 
437  *cloned_vector = *this;
438 
439  return UniquePtr<NumericVector<T> >(cloned_vector);
440 }
const Parallel::Communicator & comm() const
template<typename T >
void libMesh::LaspackVector< T >::close ( )
inlinevirtual

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

Implements libMesh::NumericVector< T >.

Definition at line 377 of file laspack_vector.h.

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

378 {
379  libmesh_assert (this->initialized());
380 
381 #ifndef NDEBUG
382  this->_is_closed = true;
383 #endif
384 }
virtual bool initialized() const
libmesh_assert(j)
const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
A reference to the Parallel::Communicator object used by this mesh.

Definition at line 87 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_petsc_snes_jacobian(), libMesh::__libmesh_petsc_snes_postcheck(), libMesh::__libmesh_petsc_snes_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::ParmetisPartitioner::_do_repartition(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::EigenSparseLinearSolver< T >::adjoint_solve(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::DofMap::attach_matrix(), libMesh::Parallel::BinSorter< KeyType, IdxType >::binsort(), libMesh::Parallel::Sort< KeyType, IdxType >::binsort(), libMesh::MeshCommunication::broadcast(), libMesh::SparseMatrix< T >::build(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::Parallel::Histogram< KeyType, IdxType >::build_histogram(), libMesh::PetscNonlinearSolver< T >::build_mat_null_space(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::DistributedVector< T >::clone(), libMesh::EigenSparseVector< T >::clone(), libMesh::LaspackVector< T >::clone(), libMesh::EpetraVector< T >::clone(), libMesh::PetscVector< T >::clone(), libMesh::EpetraVector< T >::close(), libMesh::Parallel::Sort< KeyType, IdxType >::communicate_bins(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::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::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::EquationSystems::get_solution(), libMesh::LocationMap< T >::init(), libMesh::PetscDiffSolver::init(), libMesh::TimeSolver::init(), libMesh::TopologyMap::init(), libMesh::TaoOptimizationSolver< T >::init(), libMesh::PetscNonlinearSolver< T >::init(), libMesh::DistributedVector< T >::init(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::ParmetisPartitioner::initialize(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::DistributedVector< T >::max(), libMesh::FEMSystem::mesh_position_set(), libMesh::MeshSerializer::MeshSerializer(), libMesh::DistributedVector< T >::min(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MeshTools::paranoid_n_levels(), libMesh::Partitioner::partition(), libMesh::MetisPartitioner::partition_range(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::SparseMatrix< T >::print(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Partitioner::set_parent_processor_ids(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::Parallel::Sort< KeyType, IdxType >::sort(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::Parallel::sync_element_data_by_parent_id(), libMesh::Parallel::sync_node_data_by_element_id(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::DistributedVector< T >::zero_clone(), libMesh::LaspackVector< T >::zero_clone(), libMesh::EigenSparseVector< T >::zero_clone(), libMesh::EpetraVector< T >::zero_clone(), and libMesh::PetscVector< T >::zero_clone().

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

Definition at line 116 of file numeric_vector.C.

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

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

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

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)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
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 685 of file numeric_vector.h.

687  {
688  libmesh_not_implemented();
689  }
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited
template<typename T >
T libMesh::LaspackVector< T >::dot ( const NumericVector< T > &  v) const
virtual
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, libMesh::initialized(), and libMesh::libmesh_assert().

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
libmesh_assert(j)
template<typename T>
virtual T libMesh::NumericVector< T >::el ( const numeric_index_type  i) const
inlinevirtualinherited
Returns
(*this)(i).

Definition at line 339 of file numeric_vector.h.

339 { return (*this)(i); }
void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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

Definition at line 101 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

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

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

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

469 {
470  libmesh_assert (this->initialized());
471 
472  return 0;
473 }
virtual bool initialized() const
libmesh_assert(j)
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 808 of file numeric_vector.h.

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

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

824 {
825  const std::size_t num = index.size();
826  values.resize(num);
827  if (!num)
828  return;
829 
830  this->get(index, &values[0]);
831 }
std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

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

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
template<typename T>
int libMesh::NumericVector< T >::global_relative_compare ( const NumericVector< T > &  other_vector,
const Real  threshold = TOLERANCE 
) const
virtualinherited
Returns
-1 when this is equivalent to other_vector (up to the given 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 181 of file numeric_vector.C.

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

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

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

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

Definition at line 185 of file reference_counter.h.

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

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

186 {
187  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
188  std::pair<unsigned int, unsigned int> & p = _counts[name];
189 
190  p.first++;
191 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

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

Definition at line 198 of file reference_counter.h.

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

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

199 {
200  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
201  std::pair<unsigned int, unsigned int> & p = _counts[name];
202 
203  p.second++;
204 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
template<typename T>
virtual void libMesh::LaspackVector< T >::init ( const numeric_index_type  ,
const numeric_index_type  ,
const bool  = false,
const ParallelType  = AUTOMATIC 
)
virtual

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 >::init(), libMesh::LaspackVector< T >::LaspackVector(), and libMesh::LaspackVector< T >::~LaspackVector().

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

Call init() with n_local = N.

Implements libMesh::NumericVector< T >.

Definition at line 342 of file laspack_vector.h.

References libMesh::LaspackVector< T >::init(), libMesh::libmesh_assert(), and libMesh::libmesh_dbg_var().

345 {
346  this->init(n,n,fast,ptype);
347 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
template<typename T>
virtual void libMesh::LaspackVector< T >::init ( const numeric_index_type  ,
const numeric_index_type  ,
const std::vector< numeric_index_type > &  ,
const bool  = false,
const ParallelType  = AUTOMATIC 
)
virtual

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

Implements libMesh::NumericVector< T >.

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

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

Implements libMesh::NumericVector< T >.

Definition at line 367 of file laspack_vector.h.

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

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

Definition at line 135 of file numeric_vector.h.

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

135 { return _is_initialized; }
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 94 of file numeric_vector.C.

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

96 {
97  for (numeric_index_type i=0; i<dof_indices.size(); i++)
98  this->set (dof_indices[i], v[i]);
99 }
dof_id_type numeric_index_type
Definition: id_types.h:92
template<typename T>
void libMesh::NumericVector< T >::insert ( const std::vector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

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

Definition at line 861 of file numeric_vector.h.

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

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

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

Definition at line 104 of file numeric_vector.C.

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

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

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

Definition at line 873 of file numeric_vector.h.

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

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

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

Definition at line 885 of file numeric_vector.h.

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

887 {
888  libmesh_assert(v.size() == dof_indices.size());
889  if (!v.empty())
890  this->insert(&v(0), dof_indices);
891 }
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)
libmesh_assert(j)
template<typename T >
Real libMesh::LaspackVector< T >::l1_norm ( ) const
virtual
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(), libMesh::libmesh_assert(), 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
libmesh_assert(j)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template<typename T >
Real libMesh::LaspackVector< T >::l2_norm ( ) const
virtual
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(), libMesh::libmesh_assert(), 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
libmesh_assert(j)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template<typename T >
numeric_index_type libMesh::LaspackVector< T >::last_local_index ( ) const
inlinevirtual
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 479 of file laspack_vector.h.

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

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

480 {
481  libmesh_assert (this->initialized());
482 
483  return this->size();
484 }
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type size() const libmesh_override
template<typename T >
Real libMesh::LaspackVector< T >::linfty_norm ( ) const
virtual
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(), libMesh::libmesh_assert(), 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
libmesh_assert(j)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 148 of file numeric_vector.C.

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

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

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

Implements libMesh::NumericVector< T >.

Definition at line 457 of file laspack_vector.h.

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

458 {
459  libmesh_assert (this->initialized());
460 
461  return this->size();
462 }
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::localize ( std::vector< T > &  v_local) const
virtual

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

Implements libMesh::NumericVector< T >.

Definition at line 395 of file laspack_vector.C.

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

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 }
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::localize ( NumericVector< T > &  v_local) const
virtual

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

Implements libMesh::NumericVector< T >.

Definition at line 335 of file laspack_vector.C.

References libMesh::libmesh_assert(), libMesh::libmesh_dbg_var(), libMesh::LaspackVector< T >::localize(), and libMesh::LaspackVector< T >::size().

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 }
libmesh_assert(j)
template<typename T>
virtual void libMesh::LaspackVector< T >::localize ( NumericVector< T > &  v_local,
const std::vector< numeric_index_type > &  send_list 
) const
virtual

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

Implements libMesh::NumericVector< T >.

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

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

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.

References libMesh::libmesh_dbg_var(), and libMesh::LaspackVector< T >::localize().

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
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 
)
virtual

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

Implements libMesh::NumericVector< T >.

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

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

Implements libMesh::NumericVector< T >.

Definition at line 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 libmesh_override
template<typename T >
Real libMesh::LaspackVector< T >::max ( ) const
virtual
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_assert(), 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 bool initialized() const
long double max(long double a, double b)
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template<typename T >
Real libMesh::LaspackVector< T >::min ( ) const
virtual
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_assert(), 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 bool initialized() const
long double max(long double a, double b)
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
long double min(long double a, double b)
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited
processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
The number of processors in the group.

Definition at line 93 of file parallel_object.h.

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

Referenced by libMesh::ParmetisPartitioner::_do_repartition(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::MeshCommunication::broadcast(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), libMesh::MeshBase::get_info(), libMesh::EquationSystems::get_solution(), libMesh::DistributedVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::SparseMatrix< T >::print(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::CheckpointIO::read_nodes(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::Partitioner::repartition(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Parallel::Sort< KeyType, IdxType >::sort(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

94  { return cast_int<processor_id_type>(_communicator.size()); }
unsigned int size() const
Definition: parallel.h:722
const Parallel::Communicator & _communicator
template<typename T >
T libMesh::LaspackVector< T >::operator() ( const numeric_index_type  i) const
inlinevirtual
Returns
A copy of the ith entry of the vector.

Implements libMesh::NumericVector< T >.

Definition at line 522 of file laspack_vector.h.

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

523 {
524  libmesh_assert (this->initialized());
525  libmesh_assert ( ((i >= this->first_local_index()) &&
526  (i < this->last_local_index())) );
527 
528 
529  return static_cast<T>(V_GetCmp(const_cast<QVector*>(&_vec), i+1));
530 }
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type last_local_index() const libmesh_override
virtual numeric_index_type first_local_index() const libmesh_override
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 384 of file numeric_vector.h.

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

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(), and libMesh::libmesh_assert().

85 {
86  libmesh_assert (this->closed());
87 
88  this->add(1., v);
89 
90  return *this;
91 }
virtual bool closed() const
libmesh_assert(j)
virtual void add(const numeric_index_type i, const T value) libmesh_override
template<typename T >
NumericVector< T > & libMesh::LaspackVector< T >::operator-= ( const NumericVector< T > &  v)
virtual

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(), and libMesh::libmesh_assert().

98 {
99  libmesh_assert (this->closed());
100 
101  this->add(-1., v);
102 
103  return *this;
104 }
virtual bool closed() const
libmesh_assert(j)
virtual void add(const numeric_index_type i, const T value) libmesh_override
template<typename T >
NumericVector< T > & libMesh::LaspackVector< T >::operator/= ( NumericVector< T > &  )
virtual

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 }
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
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 393 of file numeric_vector.h.

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

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

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

265 {
266  libmesh_assert (this->initialized());
267  libmesh_assert (this->closed());
268 
269  V_SetAllCmp (&_vec, s);
270 
271  return *this;
272 }
virtual bool closed() const
virtual bool initialized() const
libmesh_assert(j)
template<typename T >
NumericVector< T > & libMesh::LaspackVector< T >::operator= ( const NumericVector< T > &  v)
virtual

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

Returns
A reference to *this as the base type.

Implements libMesh::NumericVector< T >.

Definition at line 278 of file laspack_vector.C.

References libMesh::libmesh_assert().

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 }
libmesh_assert(j)
template<typename T >
LaspackVector< T > & libMesh::LaspackVector< T >::operator= ( const LaspackVector< T > &  v)

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

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(), libMesh::libmesh_assert(), 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 bool initialized() const
libmesh_assert(j)
virtual numeric_index_type size() const libmesh_override
template<typename T >
NumericVector< T > & libMesh::LaspackVector< T >::operator= ( const std::vector< T > &  v)
virtual

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 }
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::pointwise_mult ( const NumericVector< T > &  vec1,
const NumericVector< T > &  vec2 
)
virtual

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 }
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 918 of file numeric_vector.h.

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

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

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

Definition at line 900 of file numeric_vector.h.

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

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

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

Definition at line 955 of file numeric_vector.h.

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

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

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

Definition at line 933 of file numeric_vector.h.

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

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

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

Definition at line 88 of file reference_counter.C.

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

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

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

Print the contents of the 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 674 of file numeric_vector.h.

675  {
676  libmesh_not_implemented();
677  }
processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
The rank of this processor in the group.

Definition at line 99 of file parallel_object.h.

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

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::EquationSystems::_read_impl(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::MeshRefinement::add_node(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::MeshCommunication::broadcast(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::ParmetisPartitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::DofMap::build_sparsity(), libMesh::DistributedMesh::clear(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::DistributedMesh::delete_node(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::EquationSystems::get_solution(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::DistributedVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MetisPartitioner::partition_range(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::SparseMatrix< T >::print(), libMesh::NumericVector< T >::print_global(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::CheckpointIO::read_header(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::Parallel::Packing< Node * >::unpack(), libMesh::Parallel::Packing< Elem * >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::MeshTools::weight(), libMesh::NameBasedIO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), and libMesh::ExodusII_IO_Helper::write_timestep().

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

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 }
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::scale ( const T  factor)
virtual

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(), and libMesh::libmesh_assert().

230 {
231  libmesh_assert (this->initialized());
232 
233  Asgn_VV(&_vec, Mul_SV (factor, &_vec));
234 }
virtual bool initialized() const
libmesh_assert(j)
template<typename T >
void libMesh::LaspackVector< T >::set ( const numeric_index_type  i,
const T  value 
)
inlinevirtual

Sets v(i) = value.

Implements libMesh::NumericVector< T >.

Definition at line 490 of file laspack_vector.h.

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

491 {
492  libmesh_assert (this->initialized());
493  libmesh_assert_less (i, this->size());
494 
495  V_SetCmp (&_vec, i+1, value);
496 
497 #ifndef NDEBUG
498  this->_is_closed = false;
499 #endif
500 }
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type size() const libmesh_override
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 323 of file numeric_vector.C.

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

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

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

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

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

342 {
343  const NumericVector<T> & v = *this;
344 
345  std::set<numeric_index_type>::const_iterator it = indices.begin();
346  const std::set<numeric_index_type>::const_iterator it_end = indices.end();
347 
348  Real norm = 0;
349 
350  for (; it!=it_end; ++it)
351  norm += TensorTools::norm_sq(v(*it));
352 
353  this->comm().sum(norm);
354 
355  return std::sqrt(norm);
356 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
template<class T >
Real libMesh::NumericVector< T >::subset_linfty_norm ( const std::set< numeric_index_type > &  indices) const
virtualinherited
Returns
The maximum absolute value of the specified entries of this vector, which is the $ \ell_{\infty} $-norm of a vector.
Note
The indices must necessarily live on this processor.

Definition at line 359 of file numeric_vector.C.

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

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

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

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 bool closed() const
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const libmesh_override
template<typename T >
void libMesh::LaspackVector< T >::swap ( NumericVector< T > &  v)
inlinevirtual

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 536 of file laspack_vector.h.

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

537 {
538  LaspackVector<T> & v = cast_ref<LaspackVector<T> &>(other);
539 
540  // This is all grossly dependent on Laspack version...
541 
542  std::swap(_vec.Name, v._vec.Name);
543  std::swap(_vec.Dim, v._vec.Dim);
544  std::swap(_vec.Instance, v._vec.Instance);
545  std::swap(_vec.LockLevel, v._vec.LockLevel);
546  std::swap(_vec.Multipl, v._vec.Multipl);
547  std::swap(_vec.OwnData, v._vec.OwnData);
548 
549  // This should still be O(1), since _vec.Cmp is just a pointer to
550  // data on the heap
551 
552  std::swap(_vec.Cmp, v._vec.Cmp);
553 }
void swap(Iterator &lhs, Iterator &rhs)
template<typename T>
ParallelType& libMesh::NumericVector< T >::type ( )
inlineinherited
Returns
The type (SERIAL, PARALLEL, GHOSTED) of the vector.

Definition at line 145 of file numeric_vector.h.

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

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

Implements libMesh::NumericVector< T >.

Definition at line 406 of file laspack_vector.h.

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

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

407 {
408  libmesh_assert (this->initialized());
409  libmesh_assert (this->closed());
410 
411  V_SetAllCmp (&_vec, 0.);
412 }
virtual bool closed() const
virtual bool initialized() const
libmesh_assert(j)
template<typename T >
UniquePtr< NumericVector< T > > libMesh::LaspackVector< T >::zero_clone ( ) const
inlinevirtual
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 418 of file laspack_vector.h.

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

419 {
420  NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
421 
422  cloned_vector->init(*this);
423 
424  return UniquePtr<NumericVector<T> >(cloned_vector);
425 }
const Parallel::Communicator & comm() const

Friends And Related Function Documentation

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

Make other Laspack datatypes friends

Definition at line 235 of file laspack_vector.h.

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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

Definition at line 143 of file reference_counter.h.

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

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 137 of file reference_counter.h.

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

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

Definition at line 132 of file reference_counter.h.

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


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