Poster of Linux kernelThe best gift for a Linux geek
PetscVector

PetscVector

Section: C Library Functions (3) Updated: Thu Apr 7 2011
Local index Up
 

NAME

PetscVector -  

SYNOPSIS


#include <petsc_vector.h>

Inherits NumericVector< T >.  

Public Member Functions


PetscVector (const ParallelType type=AUTOMATIC)

PetscVector (const unsigned int n, const ParallelType type=AUTOMATIC)

PetscVector (const unsigned int n, const unsigned int n_local, const ParallelType type=AUTOMATIC)

PetscVector (const unsigned int N, const unsigned int n_local, const std::vector< unsigned int > &ghost, const ParallelType type=AUTOMATIC)

PetscVector (Vec v)

~PetscVector ()

void close ()

void clear ()

void zero ()

AutoPtr< NumericVector< T > > clone () const

void init (const unsigned int N, const unsigned int n_local, const bool fast=false, const ParallelType type=AUTOMATIC)

void init (const unsigned int N, const bool fast=false, const ParallelType type=AUTOMATIC)

virtual void init (const unsigned int, const unsigned int, const std::vector< unsigned int > &, const bool=false, const ParallelType=AUTOMATIC)

virtual void init (const NumericVector< T > &other, const bool fast=false)

NumericVector< T > & operator= (const T s)

NumericVector< T > & operator= (const NumericVector< T > &V)

PetscVector< T > & operator= (const PetscVector< T > &V)

NumericVector< T > & operator= (const std::vector< T > &v)

Real min () const

Real max () const

T sum () const

Real l1_norm () const

Real l2_norm () const

Real linfty_norm () const

unsigned int size () const

unsigned int local_size () const

unsigned int first_local_index () const

unsigned int last_local_index () const

unsigned int map_global_to_local_index (const unsigned int i) const

T operator() (const unsigned int i) const

virtual void get (const std::vector< unsigned int > &index, std::vector< T > &values) const

NumericVector< T > & operator+= (const NumericVector< T > &V)

NumericVector< T > & operator-= (const NumericVector< T > &V)

void set (const unsigned int i, const T value)

void add (const unsigned int i, const T value)

void add (const T s)

void add (const NumericVector< T > &V)

void add (const T a, const NumericVector< T > &v)

void add_vector (const std::vector< T > &v, const std::vector< unsigned int > &dof_indices)

void add_vector (const NumericVector< T > &V, const std::vector< unsigned int > &dof_indices)

void add_vector (const NumericVector< T > &V, const SparseMatrix< T > &A)

void add_vector (const DenseVector< T > &V, const std::vector< unsigned int > &dof_indices)

virtual void insert (const std::vector< T > &v, const std::vector< unsigned int > &dof_indices)

virtual void insert (const NumericVector< T > &V, const std::vector< unsigned int > &dof_indices)

virtual void insert (const DenseVector< T > &V, const std::vector< unsigned int > &dof_indices)

virtual void insert (const DenseSubVector< T > &V, const std::vector< unsigned int > &dof_indices)

void scale (const T factor)

virtual void abs ()

virtual T dot (const NumericVector< T > &V) const

void localize (std::vector< T > &v_local) const

void localize (NumericVector< T > &v_local) const

void localize (NumericVector< T > &v_local, const std::vector< unsigned int > &send_list) const

void localize (const unsigned int first_local_idx, const unsigned int last_local_idx, const std::vector< unsigned int > &send_list)

void localize_to_one (std::vector< T > &v_local, const unsigned int proc_id=0) const

virtual void pointwise_mult (const NumericVector< T > &vec1, const NumericVector< T > &vec2)

void print_matlab (const std::string name='NULL') const

virtual void create_subvector (NumericVector< T > &subvector, const std::vector< unsigned int > &rows) const

virtual void swap (NumericVector< T > &v)

Vec vec ()

template<> void localize_to_one (std::vector< Real > &v_local, const unsigned int pid) const

template<> void localize_to_one (std::vector< Complex > &v_local, const unsigned int pid) const

virtual bool initialized () const

ParallelType type () const

ParallelType & type ()

virtual bool closed () const

virtual Real subset_l1_norm (const std::set< unsigned int > &indices)

virtual Real subset_l2_norm (const std::set< unsigned int > &indices)

virtual Real subset_linfty_norm (const std::set< unsigned int > &indices)

virtual T el (const unsigned int i) const

NumericVector< T > & operator*= (const T a)

NumericVector< T > & operator/= (const T a)

void add_vector (const NumericVector< T > &v, const ShellMatrix< T > &a)

virtual int compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const

template<> int compare (const NumericVector< float > &other_vector, const Real threshold) const

template<> int compare (const NumericVector< double > &other_vector, const Real threshold) const

template<> int compare (const NumericVector< long double > &other_vector, const Real threshold) const

template<> int compare (const NumericVector< Complex > &other_vector, const Real threshold) const

virtual void print (std::ostream &os=std::cout) const

template<> void print (std::ostream &os) const

virtual void print_global (std::ostream &os=std::cout) const

template<> void print_global (std::ostream &os) const
 

Static Public Member Functions


static AutoPtr< NumericVector< T > > build (const SolverPackage solver_package=libMesh::default_solver_package())

static std::string get_info ()

static void print_info ()

static unsigned int n_objects ()
 

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
 

Static Protected Attributes


static Counts _counts

static Threads::atomic< unsigned int > _n_objects

static Threads::spin_mutex _mutex
 

Private Types


typedef std::map< unsigned int, unsigned int > GlobalToLocalMap
 

Private Member Functions


void _get_array (void) const

void _restore_array (void) const
 

Private Attributes


Vec _vec

bool _array_is_present

unsigned int _local_size

Vec _local_form

PetscScalar * _values

GlobalToLocalMap _global_to_local_map

bool _destroy_vec_on_exit
 

Friends


std::ostream & operator<< (std::ostream &os, const NumericVector< T > &v)
 

Detailed Description

 

template<typename T> class PetscVector< T >

Petsc vector. Provides a nice interface to the Petsc C-based data structures for parallel vectors.

Author:

Benjamin S. Kirk, 2002

Definition at line 59 of file petsc_vector.h.  

Member Typedef Documentation

 

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited]Data structure to log the information. The log is identified by the class name.

Definition at line 105 of file reference_counter.h.  

template<typename T> typedef std::map<unsigned int,unsigned int> PetscVector< T >::GlobalToLocalMap [private]Type for map that maps global to local ghost cells.

Definition at line 542 of file petsc_vector.h.  

Constructor & Destructor Documentation

 

template<typename T > PetscVector< T >::PetscVector (const ParallelTypetype = AUTOMATIC) [inline, explicit]Dummy-Constructor. Dimension=0

Definition at line 564 of file petsc_vector.h.

References NumericVector< T >::_type, and NumericVector< T >::type().

  : _array_is_present(false),
    _local_form(NULL),
    _values(NULL),
    _global_to_local_map(),
    _destroy_vec_on_exit(true)
{
  this->_type = type;
}
 

template<typename T > PetscVector< T >::PetscVector (const unsigned intn, const ParallelTypetype = AUTOMATIC) [inline, explicit]Constructor. Set dimension to n and initialize all elements with zero.

Definition at line 578 of file petsc_vector.h.

References PetscVector< T >::init().

  : _array_is_present(false),
    _local_form(NULL),
    _values(NULL),
    _global_to_local_map(),
    _destroy_vec_on_exit(true)
{
  this->init(n, n, false, type);
}
 

template<typename T > PetscVector< T >::PetscVector (const unsigned intn, const unsigned intn_local, const ParallelTypetype = AUTOMATIC) [inline]Constructor. Set local dimension to n_local, the global dimension to n, and initialize all elements with zero.

Definition at line 593 of file petsc_vector.h.

References PetscVector< T >::init().

  : _array_is_present(false),
    _local_form(NULL),
    _values(NULL),
    _global_to_local_map(),
    _destroy_vec_on_exit(true)
{
  this->init(n, n_local, false, type);
}
 

template<typename T > PetscVector< T >::PetscVector (const unsigned intN, const unsigned intn_local, const std::vector< unsigned int > &ghost, const ParallelTypetype = 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 609 of file petsc_vector.h.

References PetscVector< T >::init().

  : _array_is_present(false),
    _local_form(NULL),
    _values(NULL),
    _global_to_local_map(),
    _destroy_vec_on_exit(true)
{
  this->init(n, n_local, ghost, false, type);
}
 

template<typename T > PetscVector< T >::PetscVector (Vecv) [inline]Constructor. Creates a PetscVector assuming you already have a valid PETSc Vec object. In this case, v is NOT destroyed by the PetscVector constructor when this object goes out of scope. This allows ownership of v to remain with the original creator, and to simply provide additional functionality with the PetscVector.

Definition at line 628 of file petsc_vector.h.

References PetscVector< T >::_global_to_local_map, NumericVector< T >::_is_closed, NumericVector< T >::_is_initialized, NumericVector< T >::_type, PetscVector< T >::_vec, libMesh::COMM_WORLD, libMeshEnums::GHOSTED, PetscVector< T >::local_size(), libMeshEnums::PARALLEL, and libMeshEnums::SERIAL.

  : _array_is_present(false),
    _local_form(NULL),
    _values(NULL),
    _global_to_local_map(),
    _destroy_vec_on_exit(false)
{
  this->_vec = v;
  this->_is_closed = true;
  this->_is_initialized = true;

  /* We need to ask PETSc about the (local to global) ghost value
     mapping and create the inverse mapping out of it.  */
  int ierr=0, petsc_size=0, petsc_local_size=0;
  ierr = VecGetSize(_vec, &petsc_size);
  CHKERRABORT(libMesh::COMM_WORLD,ierr);
  ierr = VecGetLocalSize(_vec, &petsc_local_size);
  CHKERRABORT(libMesh::COMM_WORLD,ierr);

  /* 
petsc_local_size is the number of non-ghost values. If it equals the global size, then we are a serial vector, and there are no ghost values. */ if(petsc_size!=petsc_local_size) { ISLocalToGlobalMapping mapping = _vec->mapping; // If is a sparsely stored vector, set up our new mapping if (mapping) { const unsigned int local_size = static_cast<unsigned int>(petsc_local_size); const unsigned int ghost_begin = static_cast<unsigned int>(petsc_local_size); const unsigned int ghost_end = static_cast<unsigned int>(mapping->n); for(unsigned int i=ghost_begin; i<ghost_end; i++) _global_to_local_map[mapping->indices[i]] = i-local_size; this->_type = GHOSTED; } else this->_type = PARALLEL; } else this->_type = SERIAL; }
 

template<typename T > PetscVector< T >::~PetscVector () [inline]Destructor, deallocates memory. Made virtual to allow for derived classes to behave properly.

Definition at line 677 of file petsc_vector.h.

{
  this->clear ();
}
 

Member Function Documentation

 

template<typename T > void PetscVector< T >::_get_array (void) const [inline, private]Queries the array (and the local form if the vector is ghosted) from Petsc.

Definition at line 1146 of file petsc_vector.h.

References libMesh::COMM_WORLD, libMeshEnums::GHOSTED, and libMesh::initialized().

{
  libmesh_assert (this->initialized());
  if(!_array_is_present)
    {
      int ierr=0;
      if(this->type() != GHOSTED)
        {
          ierr = VecGetArray(_vec, &_values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
        }
      else
        {
          ierr = VecGhostGetLocalForm (_vec,&_local_form);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          ierr = VecGetArray(_local_form, &_values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
#ifndef NDEBUG
          int local_size = 0;
          ierr = VecGetLocalSize(_local_form, &local_size);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          _local_size = static_cast<unsigned int>(local_size);
#endif
        }
      _array_is_present = true;
    }
}
 

template<typename T > void PetscVector< T >::_restore_array (void) const [inline, private]Restores the array (and the local form if the vector is ghosted) to Petsc.

Definition at line 1178 of file petsc_vector.h.

References libMesh::COMM_WORLD, libMeshEnums::GHOSTED, and libMesh::initialized().

Referenced by PetscVector< T >::add(), PetscVector< T >::create_subvector(), PetscVector< T >::init(), and PetscVector< T >::operator=().

{
  libmesh_assert (this->initialized());
  if(_array_is_present)
    {
      int ierr=0;
      if(this->type() != GHOSTED)
        {
          ierr = VecRestoreArray (_vec, &_values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          _values = NULL;
        }
      else
        {
          ierr = VecRestoreArray (_local_form, &_values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          _values = NULL;
          ierr = VecGhostRestoreLocalForm (_vec,&_local_form);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          _local_form = NULL;
#ifndef NDEBUG
          _local_size = 0;
#endif
        }
      _array_is_present = false;
    }
}
 

template<typename T > void PetscVector< T >::abs () [virtual]v = abs(v)... that is, each entry in v is replaced by its absolute value.

Implements NumericVector< T >.

Definition at line 458 of file petsc_vector.C.

References libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

{
  this->_restore_array();

  int ierr = 0;

  if(this->type() != GHOSTED)
    {
      ierr = VecAbs(_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);  
    }
  else
    {
      Vec loc_vec;
      ierr = VecGhostGetLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);

      ierr = VecAbs(loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);  

      ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }
}
 

template<typename T > void PetscVector< T >::add (const unsigned inti, const Tvalue) [virtual]v(i) += value

Implements NumericVector< T >.

Definition at line 166 of file petsc_vector.C.

References libMesh::COMM_WORLD.

{
  this->_restore_array();
  libmesh_assert(i<size());
  
  int ierr=0;
  int i_val = static_cast<int>(i);
  PetscScalar petsc_value = static_cast<PetscScalar>(value);

  ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, ADD_VALUES);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  this->_is_closed = false;
}
 

template<typename T > void PetscVector< T >::add (const Ta, const NumericVector< T > &v) [virtual]$ U+=a*V $ . Simple vector addition, equal to the operator +=.

Implements NumericVector< T >.

Definition at line 315 of file petsc_vector.C.

References PetscVector< T >::_restore_array(), PetscVector< T >::_vec, libMesh::COMM_WORLD, libMeshEnums::GHOSTED, and PetscVector< T >::size().

{
  this->_restore_array();

  int ierr = 0;
  PetscScalar a = static_cast<PetscScalar>(a_in);

  // Make sure the NumericVector passed in is really a PetscVector
  const PetscVector<T>* v = libmesh_cast_ptr<const PetscVector<T>*>(&v_in);
  v->_restore_array();

  libmesh_assert(this->size() == v->size());
  
  if(this->type() != GHOSTED)
    {
#if PETSC_VERSION_LESS_THAN(2,3,0)
      // 2.2.x & earlier style
      ierr = VecAXPY(&a, v->_vec, _vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
      // 2.3.x & later style
      ierr = VecAXPY(_vec, a, v->_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif
    }
  else
    {
      Vec loc_vec;
      Vec v_loc_vec;
      ierr = VecGhostGetLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
      ierr = VecGhostGetLocalForm (v->_vec,&v_loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);

#if PETSC_VERSION_LESS_THAN(2,3,0)
      // 2.2.x & earlier style
      ierr = VecAXPY(&a, v_loc_vec, loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
      // 2.3.x & later style
      ierr = VecAXPY(loc_vec, a, v_loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif

      ierr = VecGhostRestoreLocalForm (v->_vec,&v_loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
      ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }
}
 

template<typename T > void PetscVector< T >::add (const Ts) [virtual]$U(0-LIBMESH_DIM)+=s$. Addition of s to all components. Note that s is a scalar and not a vector.

Implements NumericVector< T >.

Definition at line 242 of file petsc_vector.C.

References libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

{
  this->_restore_array();

  int ierr=0;
  PetscScalar* values;
  const PetscScalar v = static_cast<PetscScalar>(v_in);  

  if(this->type() != GHOSTED)
    {
      const int n   = static_cast<int>(this->local_size());
      const int fli = static_cast<int>(this->first_local_index());
      
      for (int i=0; i<n; i++)
        {
          ierr = VecGetArray (_vec, &values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          
          int ig = fli + i;      
          
          PetscScalar value = (values[i] + v);
          
          ierr = VecRestoreArray (_vec, &values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          
          ierr = VecSetValues (_vec, 1, &ig, &value, INSERT_VALUES);
          CHKERRABORT(libMesh::COMM_WORLD,ierr); 
        }
    }
  else
    {
      /* Vectors that include ghost values require a special
         handling.  */
      Vec loc_vec;
      ierr = VecGhostGetLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);

      int n=0;
      ierr = VecGetSize(loc_vec, &n);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
      
      for (int i=0; i<n; i++)
        {
          ierr = VecGetArray (loc_vec, &values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          
          PetscScalar value = (values[i] + v);
          
          ierr = VecRestoreArray (loc_vec, &values);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          
          ierr = VecSetValues (loc_vec, 1, &i, &value, INSERT_VALUES);
          CHKERRABORT(libMesh::COMM_WORLD,ierr); 
        }

      ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  this->_is_closed = false;
}
 

template<typename T > void PetscVector< T >::add (const NumericVector< T > &V) [virtual]$ U+=V $ . Simple vector addition, equal to the operator +=.

Implements NumericVector< T >.

Definition at line 307 of file petsc_vector.C.

{
  this->add (1., v);
}
 

template<typename T > void PetscVector< T >::add_vector (const std::vector< T > &v, const std::vector< unsigned int > &dof_indices) [virtual]$ U+=v $ where v is a std::vector<T> and you want to specify WHERE to add it

Implements NumericVector< T >.

Definition at line 184 of file petsc_vector.C.

{
  this->_restore_array();
  libmesh_assert (v.size() == dof_indices.size());

  for (unsigned int i=0; i<v.size(); i++)
    this->add (dof_indices[i], v[i]);
}
 

template<typename T > void PetscVector< T >::add_vector (const NumericVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$ U+=V $ where U and V are type NumericVector<T> and you want to specify WHERE to add the NumericVector<T> V

Implements NumericVector< T >.

Definition at line 197 of file petsc_vector.C.

References NumericVector< T >::size().

{
  libmesh_assert (V.size() == dof_indices.size());

  for (unsigned int i=0; i<V.size(); i++)
    this->add (dof_indices[i], V(i));
}
 

template<typename T > void PetscVector< T >::add_vector (const NumericVector< T > &V, const SparseMatrix< T > &A) [virtual]$U+=A*V$, add the product of a SparseMatrix A and a NumericVector V to this, where this=U.

Implements NumericVector< T >.

Definition at line 209 of file petsc_vector.C.

References PetscVector< T >::_vec, PetscMatrix< T >::close(), and libMesh::COMM_WORLD.

{
  this->_restore_array();
  // Make sure the data passed in are really of Petsc types
  const PetscVector<T>* V = libmesh_cast_ptr<const PetscVector<T>*>(&V_in);
  const PetscMatrix<T>* A = libmesh_cast_ptr<const PetscMatrix<T>*>(&A_in);

  int ierr=0;

  A->close();

  // The const_cast<> is not elegant, but it is required since PETSc
  // is not const-correct.  
  ierr = MatMultAdd(const_cast<PetscMatrix<T>*>(A)->mat(), V->_vec, _vec, _vec);
         CHKERRABORT(libMesh::COMM_WORLD,ierr); 
}
 

template<typename T> void NumericVector< T >::add_vector (const NumericVector< T > &v, const ShellMatrix< T > &a) [inherited]$U+=A*V$, add the product of a ShellMatrix A and a NumericVector V to this, where this=U.

Definition at line 253 of file numeric_vector.C.

References ShellMatrix< T >::vector_mult_add().

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

template<typename T > void PetscVector< T >::add_vector (const DenseVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$U+=V $ where U and V are type DenseVector<T> and you want to specify WHERE to add the DenseVector<T> V

Implements NumericVector< T >.

Definition at line 230 of file petsc_vector.C.

References DenseVector< T >::size().

{
  libmesh_assert (V.size() == dof_indices.size());

  for (unsigned int i=0; i<V.size(); i++)
    this->add (dof_indices[i], V(i));
}
 

template<typename T > AutoPtr< NumericVector< T > > NumericVector< T >::build (const SolverPackagesolver_package = libMesh::default_solver_package()) [static, inherited]Builds a NumericVector using the linear solver package specified by solver_package

Definition at line 41 of file numeric_vector.C.

References LASPACK_SOLVERS, libMeshEnums::PETSC_SOLVERS, and TRILINOS_SOLVERS.

Referenced by ExactErrorEstimator::estimate_error().

{
  // Build the appropriate vector
  switch (solver_package)
    {


#ifdef LIBMESH_HAVE_LASPACK
    case LASPACK_SOLVERS:
      {
        AutoPtr<NumericVector<T> > ap(new LaspackVector<T>);
        return ap;
      }
#endif


#ifdef LIBMESH_HAVE_PETSC
    case PETSC_SOLVERS:
      {
        AutoPtr<NumericVector<T> > ap(new PetscVector<T>);
        return ap;
      }
#endif


#ifdef LIBMESH_HAVE_TRILINOS
    case TRILINOS_SOLVERS:
      {
        AutoPtr<NumericVector<T> > ap(new EpetraVector<T>);
        return ap;
      }
#endif


    default:
      AutoPtr<NumericVector<T> > ap(new DistributedVector<T>);
      return ap;
    }
    
  AutoPtr<NumericVector<T> > ap(NULL);
  return ap;    
}
 

template<typename T > void PetscVector< T >::clear () [inline, virtual]Returns:

the PetscVector<T> to a pristine state.

Reimplemented from NumericVector< T >.

Definition at line 875 of file petsc_vector.h.

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

{
  if (this->initialized())
    this->_restore_array();

  if ((this->initialized()) && (this->_destroy_vec_on_exit))
    {
      int ierr=0;

      ierr = VecDestroy(_vec);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  this->_is_closed = this->_is_initialized = false;

  _global_to_local_map.clear();
}
 

template<typename T > AutoPtr< NumericVector< T > > PetscVector< T >::clone () const [inline, virtual]Creates a copy of this vector and returns it in an AutoPtr.

Implements NumericVector< T >.

Definition at line 942 of file petsc_vector.h.

{
  AutoPtr<NumericVector<T> > cloned_vector (new PetscVector<T>);

  cloned_vector->init(*this, true);

  *cloned_vector = *this;

  return cloned_vector;
}
 

template<typename T > void PetscVector< T >::close () [inline, virtual]Call the assemble functions

Implements NumericVector< T >.

Definition at line 849 of file petsc_vector.h.

References libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

Referenced by SlepcEigenSolver< T >::get_eigenpair(), PetscVector< T >::localize(), PetscLinearSolver< T >::solve(), and PetscDiffSolver::solve().

{
  this->_restore_array();
  
  int ierr=0;
  
  ierr = VecAssemblyBegin(_vec);
  CHKERRABORT(libMesh::COMM_WORLD,ierr);
  ierr = VecAssemblyEnd(_vec);
  CHKERRABORT(libMesh::COMM_WORLD,ierr);

  if(this->type() == GHOSTED)
    {
      ierr = VecGhostUpdateBegin(_vec,INSERT_VALUES,SCATTER_FORWARD);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
      ierr = VecGhostUpdateEnd(_vec,INSERT_VALUES,SCATTER_FORWARD);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  this->_is_closed = true;
}
 

template<typename T> virtual bool NumericVector< T >::closed () const [inline, virtual, inherited]Returns:

true if the vector is closed and ready for computation, false otherwise.

Definition at line 125 of file numeric_vector.h.

Referenced by DofMap::enforce_constraints_exactly(), and DofMap::max_constraint_error().

{ return _is_closed; }
 

template<typename T> virtual int NumericVector< T >::compare (const NumericVector< T > &other_vector, const Realthreshold = TOLERANCE) const [virtual, inherited]Returns:

-1 when this is equivalent to other_vector, up to the given threshold. When differences occur, the return value contains the first index where the difference exceeded the threshold. When no threshold is given, the libMesh TOLERANCE is used.

 

template<> int NumericVector< float >::compare (const NumericVector< float > &other_vector, const Realthreshold) const [inherited]

Definition at line 89 of file numeric_vector.C.

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

{
  libmesh_assert (this->initialized());
  libmesh_assert (other_vector.initialized());
  libmesh_assert (this->first_local_index() == other_vector.first_local_index());
  libmesh_assert (this->last_local_index()  == other_vector.last_local_index());

  int rvalue     = -1;
  unsigned int i = first_local_index();

  do
    {
      if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
        rvalue = i;
      else
        i++;
    }
  while (rvalue==-1 && i<last_local_index());

  return rvalue;
}
 

template<> int NumericVector< double >::compare (const NumericVector< double > &other_vector, const Realthreshold) const [inherited]

Definition at line 114 of file numeric_vector.C.

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

{
  libmesh_assert (this->initialized());
  libmesh_assert (other_vector.initialized());
  libmesh_assert (this->first_local_index() == other_vector.first_local_index());
  libmesh_assert (this->last_local_index()  == other_vector.last_local_index());

  int rvalue     = -1;
  unsigned int i = first_local_index();

  do
    {
      if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
        rvalue = i;
      else
        i++;
    }
  while (rvalue==-1 && i<last_local_index());

  return rvalue;
}
 

template<> int NumericVector< long double >::compare (const NumericVector< long double > &other_vector, const Realthreshold) const [inherited]

Definition at line 140 of file numeric_vector.C.

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

{
  libmesh_assert (this->initialized());
  libmesh_assert (other_vector.initialized());
  libmesh_assert (this->first_local_index() == other_vector.first_local_index());
  libmesh_assert (this->last_local_index()  == other_vector.last_local_index());

  int rvalue     = -1;
  unsigned int i = first_local_index();

  do
    {
      if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
        rvalue = i;
      else
        i++;
    }
  while (rvalue==-1 && i<last_local_index());

  return rvalue;
}
 

template<> int NumericVector< Complex >::compare (const NumericVector< Complex > &other_vector, const Realthreshold) const [inherited]

Definition at line 167 of file numeric_vector.C.

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

{
  libmesh_assert (this->initialized());
  libmesh_assert (other_vector.initialized());
  libmesh_assert (this->first_local_index() == other_vector.first_local_index());
  libmesh_assert (this->last_local_index()  == other_vector.last_local_index());

  int rvalue     = -1;
  unsigned int i = first_local_index();

  do
    {
      if (( std::abs( (*this)(i).real() - other_vector(i).real() ) > threshold ) ||
          ( std::abs( (*this)(i).imag() - other_vector(i).imag() ) > threshold ))
        rvalue = i;
      else
        i++;
    }
  while (rvalue==-1 && i<this->last_local_index());

  return rvalue;
}
 

template<typename T > void PetscVector< T >::create_subvector (NumericVector< T > &subvector, const std::vector< unsigned int > &rows) const [virtual]Creates a 'subvector' from this vector using the rows indices of the 'rows' array.

Reimplemented from NumericVector< T >.

Definition at line 1195 of file petsc_vector.C.

References NumericVector< T >::_is_initialized, PetscVector< T >::_restore_array(), PetscVector< T >::_vec, libMesh::COMM_WORLD, MeshTools::Generation::Private::idx(), NumericVector< T >::initialized(), and Utility::iota().

{
  this->_restore_array();

  // PETSc data structures
  IS parent_is, subvector_is;
  VecScatter scatter;
  int ierr = 0;
  
  // Make sure the passed in subvector is really a PetscVector
  PetscVector<T>* petsc_subvector = libmesh_cast_ptr<PetscVector<T>*>(&subvector);
  
  // If the petsc_subvector is already initialized, we assume that the
  // user has already allocated the *correct* amount of space for it.
  // If not, we use the appropriate PETSc routines to initialize it.
  if (!petsc_subvector->initialized())
    {
      // Initialize the petsc_subvector to have enough space to hold
      // the entries which will be scattered into it.  Note: such an
      // init() function (where we let PETSc decide the number of local
      // entries) is not currently offered by the PetscVector
      // class.  Should we differentiate here between sequential and
      // parallel vector creation based on libMesh::n_processors() ?
      ierr = VecCreateMPI(libMesh::COMM_WORLD,
                          PETSC_DECIDE,          // n_local
                          rows.size(),           // n_global
                          &(petsc_subvector->_vec)); CHKERRABORT(libMesh::COMM_WORLD,ierr);

      ierr = VecSetFromOptions (petsc_subvector->_vec); CHKERRABORT(libMesh::COMM_WORLD,ierr);

      // Mark the subvector as initialized
      petsc_subvector->_is_initialized = true;
    }
  else
    {
      petsc_subvector->_restore_array();
    }
  
  // Use iota to fill an array with entries [0,1,2,3,4,...rows.size()]
  std::vector<int> idx(rows.size());
  Utility::iota (idx.begin(), idx.end(), 0);

  // Construct index sets
  ierr = ISCreateGeneral(libMesh::COMM_WORLD,
                         rows.size(),
                         (int*) &rows[0],
                         &parent_is); CHKERRABORT(libMesh::COMM_WORLD,ierr);

  ierr = ISCreateGeneral(libMesh::COMM_WORLD,
                         rows.size(),
                         (int*) &idx[0],
                         &subvector_is); CHKERRABORT(libMesh::COMM_WORLD,ierr);

  // Construct the scatter object
  ierr = VecScatterCreate(this->_vec,
                          parent_is,
                          petsc_subvector->_vec,
                          subvector_is,
                          &scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr);

  // Actually perform the scatter
#if PETSC_VERSION_LESS_THAN(2,3,3)
  ierr = VecScatterBegin(this->_vec,
                         petsc_subvector->_vec,
                         INSERT_VALUES,
                         SCATTER_FORWARD,
                         scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr);

  ierr = VecScatterEnd(this->_vec,
                       petsc_subvector->_vec,
                       INSERT_VALUES,
                       SCATTER_FORWARD,
                       scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
  // API argument order change in PETSc 2.3.3
  ierr = VecScatterBegin(scatter,
                         this->_vec,
                         petsc_subvector->_vec,
                         INSERT_VALUES,
                         SCATTER_FORWARD); CHKERRABORT(libMesh::COMM_WORLD,ierr);

  ierr = VecScatterEnd(scatter,
                       this->_vec,
                       petsc_subvector->_vec,
                       INSERT_VALUES,
                       SCATTER_FORWARD); CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif
  
  // Clean up 
  ierr = ISDestroy(parent_is);       CHKERRABORT(libMesh::COMM_WORLD,ierr);
  ierr = ISDestroy(subvector_is);    CHKERRABORT(libMesh::COMM_WORLD,ierr);
  ierr = VecScatterDestroy(scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr); 

}
 

template<typename T > T PetscVector< T >::dot (const NumericVector< T > &V) const [virtual]Computes the dot product, p = U.V

Implements NumericVector< T >.

Definition at line 484 of file petsc_vector.C.

References PetscVector< T >::_vec, and libMesh::COMM_WORLD.

{
  this->_restore_array();

  // Error flag
  int ierr = 0;
  
  // Return value
  PetscScalar value=0.;
  
  // Make sure the NumericVector passed in is really a PetscVector
  const PetscVector<T>* v = libmesh_cast_ptr<const PetscVector<T>*>(&V);

  // 2.3.x (at least) style.  Untested for previous versions.
  ierr = VecDot(this->_vec, v->_vec, &value);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  return static_cast<T>(value);
}
 

template<typename T> virtual T NumericVector< T >::el (const unsigned inti) const [inline, virtual, inherited]Returns:

the element U(i)

Definition at line 314 of file numeric_vector.h.

{ return (*this)(i); }
 

template<typename T > unsigned int PetscVector< T >::first_local_index () const [inline, virtual]Returns:

the index of the first vector element actually stored on this processor

Implements NumericVector< T >.

Definition at line 992 of file petsc_vector.h.

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

{
  libmesh_assert (this->initialized());
  
  int ierr=0, petsc_first=0, petsc_last=0;
  
  ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  return static_cast<unsigned int>(petsc_first);
}
 

template<typename T > void PetscVector< T >::get (const std::vector< unsigned int > &index, std::vector< T > &values) const [inline, virtual]Access multiple components at once. Overloaded method that should be faster (probably much faster) than calling operator() individually for each index.

Reimplemented from NumericVector< T >.

Definition at line 1068 of file petsc_vector.h.

References libMeshEnums::GHOSTED.

{
  this->_get_array();

  const unsigned int num = index.size();
  values.resize(num);

  for(unsigned int i=0; i<num; i++)
    {
      const unsigned int local_index = this->map_global_to_local_index(index[i]);
#ifndef NDEBUG
      if(this->type() == GHOSTED)
        {
          libmesh_assert(local_index<_local_size);
        }
#endif
      values[i] = static_cast<T>(_values[local_index]);
    }
}
 

std::string ReferenceCounter::get_info () [static, inherited]Gets a string containing the reference information.

Definition at line 45 of file reference_counter.C.

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

Referenced by ReferenceCounter::print_info().

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)

  std::ostringstream out;
  
  out << '
      << ' ---------------------------------------------------------------------------- 
      << '| Reference count information                                                |
      << ' ---------------------------------------------------------------------------- ;
  
  for (Counts::iterator it = _counts.begin();
       it != _counts.end(); ++it)
    {
      const std::string name(it->first);
      const unsigned int creations    = it->second.first;
      const unsigned int destructions = it->second.second;

      out << '| ' << name << ' reference count information:
          << '|  Creations:    ' << creations    << '
          << '|  Destructions: ' << destructions << ';
    }
  
  out << ' ---------------------------------------------------------------------------- ;

  return out.str();

#else

  return '';
  
#endif
}
 

void ReferenceCounter::increment_constructor_count (const std::string &name) [inline, protected, inherited]Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 149 of file reference_counter.h.

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

Referenced by ReferenceCountedObject< Value >::ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.first++;
}
 

void ReferenceCounter::increment_destructor_count (const std::string &name) [inline, protected, inherited]Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 167 of file reference_counter.h.

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

Referenced by ReferenceCountedObject< Value >::~ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.second++;
}
 

template<typename T > void PetscVector< T >::init (const NumericVector< T > &other, const boolfast = false) [inline, virtual]Creates a vector that has the same dimension and storage type as other, including ghost dofs.

Implements NumericVector< T >.

Definition at line 813 of file petsc_vector.h.

References PetscVector< T >::_global_to_local_map, NumericVector< T >::_is_closed, NumericVector< T >::_is_initialized, libMesh::libMeshPrivateData::_is_initialized, PetscVector< T >::_restore_array(), NumericVector< T >::_type, PetscVector< T >::_vec, libMesh::COMM_WORLD, NumericVector< T >::initialized(), libMesh::initialized(), PetscVector< T >::size(), and libMesh::zero.

{
  // Clear initialized vectors 
  if (this->initialized())
    this->clear();

  const PetscVector<T>& v = libmesh_cast_ref<const PetscVector<T>&>(other);

  // Other vector should restore array.
  if(v.initialized())
    {
      v._restore_array();
    }

  this->_global_to_local_map = v._global_to_local_map;
  this->_is_closed      = v._is_closed;
  this->_is_initialized = v._is_initialized;
  this->_type = v._type;

  if (v.size() != 0)
    {
      int ierr = 0;

      ierr = VecDuplicate (v._vec, &this->_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }
  
  if (fast == false)
    this->zero ();
}
 

template<typename T > void PetscVector< T >::init (const unsigned intN, const unsigned intn_local, const boolfast = false, const ParallelTypetype = AUTOMATIC) [inline, virtual]Change the dimension of the vector to N. The reserved memory for this vector remains unchanged if possible, to make things faster, but this may waste some memory, so take this in the back of your head. However, if N==0 all memory is freed, i.e. if you want to resize the vector and release the memory not needed, you have to first call init(0) and then init(N). This cited behaviour is analogous to that of the STL containers.

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

Implements NumericVector< T >.

Definition at line 686 of file petsc_vector.h.

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

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

{
  int ierr=0;
  int petsc_n=static_cast<int>(n);
  int petsc_n_local=static_cast<int>(n_local);


  // Clear initialized vectors 
  if (this->initialized())
    this->clear();

  if (type == AUTOMATIC)
    {
      if (n == n_local)
        this->_type = SERIAL;
      else
        this->_type = PARALLEL;
    }
  else
    this->_type = type;

  libmesh_assert ((this->_type==SERIAL && n==n_local) ||
                  this->_type==PARALLEL);
  
  // create a sequential vector if on only 1 processor 
  if (this->_type == SERIAL)
    {
      ierr = VecCreateSeq (PETSC_COMM_SELF, petsc_n, &_vec);
             CHKERRABORT(PETSC_COMM_SELF,ierr);
      
      ierr = VecSetFromOptions (_vec);
             CHKERRABORT(PETSC_COMM_SELF,ierr);
    }
  // otherwise create an MPI-enabled vector
  else if (this->_type == PARALLEL)
    {
      libmesh_assert (n_local <= n);
      
      ierr = VecCreateMPI (libMesh::COMM_WORLD, petsc_n_local, petsc_n,
                           &_vec);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
      
      ierr = VecSetFromOptions (_vec);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }  
  else
    libmesh_error();
  
  this->_is_initialized = true;
//  this->_is_closed = true;
  
  
  if (fast == false)
    this->zero ();
}
 

template<typename T > void PetscVector< T >::init (const unsigned intN, const boolfast = false, const ParallelTypetype = AUTOMATIC) [inline, virtual]call init with n_local = N,

Implements NumericVector< T >.

Definition at line 749 of file petsc_vector.h.

References libMesh::init().

{
  this->init(n,n,fast,type);
}
 

template<typename T > void PetscVector< T >::init (const unsigned intn, const unsigned intn_local, const std::vector< unsigned int > &ghost, const boolfast = false, const ParallelTypetype = AUTOMATIC) [inline, virtual]Create a vector that holds tha local indices plus those specified in the ghost argument.

Implements NumericVector< T >.

Definition at line 760 of file petsc_vector.h.

References libMesh::libMeshPrivateData::_is_initialized, libMeshEnums::AUTOMATIC, libMesh::COMM_WORLD, libMeshEnums::GHOSTED, libMesh::initialized(), and libMesh::zero.

{
  int ierr=0;
  int petsc_n=static_cast<int>(n);
  int petsc_n_local=static_cast<int>(n_local);
  int petsc_n_ghost=static_cast<int>(ghost.size());

  // If the mesh is disjoint, the following assertion will fail.
  // If the mesh is not disjoint, every processor will either have
  // all the dofs, none of the dofs, or some non-zero dofs at the
  // boundary between processors.
  libmesh_assert(n_local == 0 || n_local == n || !ghost.empty());

  int* petsc_ghost = ghost.empty() ? PETSC_NULL :
    const_cast<int*>(reinterpret_cast<const int*>(&ghost[0]));

  // Clear initialized vectors 
  if (this->initialized())
    this->clear();

  libmesh_assert(type == AUTOMATIC || type == GHOSTED);
  this->_type = GHOSTED;

  /* Make the global-to-local ghost cell map.  */
  for(unsigned int i=0; i<ghost.size(); i++)
    {
      _global_to_local_map[ghost[i]] = i;
    }

  /* Create vector.  */
  ierr = VecCreateGhost (libMesh::COMM_WORLD, petsc_n_local, petsc_n,
                         petsc_n_ghost, petsc_ghost,
                         &_vec);
  CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  ierr = VecSetFromOptions (_vec);
  CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  this->_is_initialized = true;
//  this->_is_closed = true;
  
  if (fast == false)
    this->zero ();
}
 

template<typename T> virtual bool NumericVector< T >::initialized () const [inline, virtual, inherited]Returns:

true if the vector has been initialized, false otherwise.

Definition at line 109 of file numeric_vector.h.

Referenced by ImplicitSystem::assemble(), ExplicitSystem::assemble_qoi(), ExplicitSystem::assemble_qoi_derivative(), NumericVector< T >::compare(), PetscVector< T >::create_subvector(), and PetscVector< T >::init().

{ return _is_initialized; }
 

template<typename T > void PetscVector< T >::insert (const std::vector< T > &v, const std::vector< unsigned int > &dof_indices) [virtual]$ U=v $ where v is a std::vector<T> and you want to specify WHERE to insert it

Implements NumericVector< T >.

Definition at line 369 of file petsc_vector.C.

{
  libmesh_assert (v.size() == dof_indices.size());

  for (unsigned int i=0; i<v.size(); i++)
    this->set (dof_indices[i], v[i]);
}
 

template<typename T > void PetscVector< T >::insert (const NumericVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$U=V$, where U and V are type NumericVector<T> and you want to specify WHERE to insert the NumericVector<T> V

Implements NumericVector< T >.

Definition at line 381 of file petsc_vector.C.

References NumericVector< T >::size().

{
  libmesh_assert (V.size() == dof_indices.size());

  for (unsigned int i=0; i<V.size(); i++)
    this->set (dof_indices[i], V(i));
}
 

template<typename T > void PetscVector< T >::insert (const DenseSubVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$ U=V $ where V is type DenseSubVector<T> and you want to specify WHERE to insert it

Implements NumericVector< T >.

Definition at line 405 of file petsc_vector.C.

References DenseVectorBase< T >::size().

{
  libmesh_assert (V.size() == dof_indices.size());

  for (unsigned int i=0; i<V.size(); i++)
    this->set (dof_indices[i], V(i));
}
 

template<typename T > void PetscVector< T >::insert (const DenseVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$ U=V $ where V is type DenseVector<T> and you want to specify WHERE to insert it

Implements NumericVector< T >.

Definition at line 393 of file petsc_vector.C.

References DenseVector< T >::size().

{
  libmesh_assert (V.size() == dof_indices.size());

  for (unsigned int i=0; i<V.size(); i++)
    this->set (dof_indices[i], V(i));
}
 

template<typename T > Real PetscVector< T >::l1_norm () const [virtual]Returns:

the $l_1$-norm of the vector, i.e. the sum of the absolute values.

Implements NumericVector< T >.

Definition at line 67 of file petsc_vector.C.

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

{
  this->_restore_array();
  libmesh_assert(this->closed());
  
  int ierr=0;
  PetscReal value=0.;
  
  ierr = VecNorm (_vec, NORM_1, &value);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  return static_cast<Real>(value);
}
 

template<typename T > Real PetscVector< T >::l2_norm () const [virtual]Returns:

the $l_2$-norm of the vector, i.e. the square root of the sum of the squares of the elements.

Implements NumericVector< T >.

Definition at line 84 of file petsc_vector.C.

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

{
  this->_restore_array();
  libmesh_assert(this->closed());
  
  int ierr=0;
  PetscReal value=0.;
  
  ierr = VecNorm (_vec, NORM_2, &value);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  return static_cast<Real>(value);
}
 

template<typename T > unsigned int PetscVector< T >::last_local_index () const [inline, virtual]Returns:

the index of the last vector element actually stored on this processor

Implements NumericVector< T >.

Definition at line 1008 of file petsc_vector.h.

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

{
  libmesh_assert (this->initialized());
  
  int ierr=0, petsc_first=0, petsc_last=0;
  
  ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  return static_cast<unsigned int>(petsc_last);
}
 

template<typename T > Real PetscVector< T >::linfty_norm () const [virtual]Returns:

the maximum absolute value of the elements of this vector, which is the $l_infty$-norm of a vector.

Implements NumericVector< T >.

Definition at line 102 of file petsc_vector.C.

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

{
  this->_restore_array();
  libmesh_assert(this->closed());
  
  int ierr=0;
  PetscReal value=0.;
  
  ierr = VecNorm (_vec, NORM_INFINITY, &value);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  return static_cast<Real>(value);
}
 

template<typename T > unsigned int PetscVector< T >::local_size () const [inline, virtual]Returns:

the local size of the vector (index_stop-index_start)

Implements NumericVector< T >.

Definition at line 976 of file petsc_vector.h.

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

Referenced by PetscVector< T >::operator=(), and PetscVector< T >::PetscVector().

{
  libmesh_assert (this->initialized());
  
  int ierr=0, petsc_size=0;
  
  ierr = VecGetLocalSize(_vec, &petsc_size);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  return static_cast<unsigned int>(petsc_size);
}
 

template<typename T > void PetscVector< T >::localize (std::vector< T > &v_local) const [virtual]Creates a copy of the global vector in the local vector v_local.

Implements NumericVector< T >.

Definition at line 901 of file petsc_vector.C.

References libMesh::COMM_WORLD.

Referenced by PetscVector< T >::localize().

{
  this->_restore_array();

  // This function must be run on all processors at once
  parallel_only();

  int ierr=0;
  const int n = this->size();
  const int nl = this->local_size();
  PetscScalar *values;

  v_local.clear();
  v_local.resize(n, 0.);

  ierr = VecGetArray (_vec, &values);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  unsigned int ioff = first_local_index();

  for (int i=0; i<nl; i++)
    v_local[i+ioff] = static_cast<T>(values[i]);

  ierr = VecRestoreArray (_vec, &values);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  Parallel::sum(v_local);
}
 

template<typename T > void PetscVector< T >::localize (NumericVector< T > &v_local) const [virtual]Same, but fills a NumericVector<T> instead of a std::vector.

Implements NumericVector< T >.

Definition at line 670 of file petsc_vector.C.

References PetscVector< T >::_vec, PetscVector< T >::close(), libMesh::COMM_WORLD, libMeshEnums::GHOSTED, MeshTools::Generation::Private::idx(), Utility::iota(), PetscVector< T >::size(), and NumericVector< T >::type().

{
  this->_restore_array();

  // Make sure the NumericVector passed in is really a PetscVector
  PetscVector<T>* v_local = libmesh_cast_ptr<PetscVector<T>*>(&v_local_in);

  libmesh_assert (v_local != NULL);
  libmesh_assert (v_local->size() == this->size());

  int ierr = 0;
  const int n = this->size();

  IS is;
  VecScatter scatter;

  // Create idx, idx[i] = i;
  std::vector<int> idx(n); Utility::iota (idx.begin(), idx.end(), 0);

  // Create the index set & scatter object
  ierr = ISCreateGeneral(libMesh::COMM_WORLD, n, &idx[0], &is);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  ierr = VecScatterCreate(_vec,          is,
                          v_local->_vec, is,
                          &scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  // Perform the scatter
#if PETSC_VERSION_LESS_THAN(2,3,3)
         
  ierr = VecScatterBegin(_vec, v_local->_vec, INSERT_VALUES,
                         SCATTER_FORWARD, scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  ierr = VecScatterEnd  (_vec, v_local->_vec, INSERT_VALUES,
                         SCATTER_FORWARD, scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
  // API argument order change in PETSc 2.3.3
  ierr = VecScatterBegin(scatter, _vec, v_local->_vec,
                         INSERT_VALUES, SCATTER_FORWARD);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  ierr = VecScatterEnd  (scatter, _vec, v_local->_vec,
                         INSERT_VALUES, SCATTER_FORWARD);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif

  // Clean up
  ierr = ISDestroy (is);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  ierr = VecScatterDestroy(scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  // Make sure ghost dofs are up to date
  if (v_local->type() == GHOSTED)
    v_local->close();
}
 

template<typename T > void PetscVector< T >::localize (NumericVector< T > &v_local, const std::vector< unsigned int > &send_list) const [virtual]Creates a local vector v_local containing only information relevant to this processor, as defined by the send_list.

Implements NumericVector< T >.

Definition at line 734 of file petsc_vector.C.

References PetscVector< T >::_vec, PetscVector< T >::close(), libMesh::COMM_WORLD, libMeshEnums::GHOSTED, MeshTools::Generation::Private::idx(), PetscVector< T >::size(), and NumericVector< T >::type().

{
  this->_restore_array();

  // Make sure the NumericVector passed in is really a PetscVector
  PetscVector<T>* v_local = libmesh_cast_ptr<PetscVector<T>*>(&v_local_in);

  libmesh_assert (v_local != NULL);
  libmesh_assert (v_local->size() == this->size());
  libmesh_assert (send_list.size()     <= v_local->size());
  
  int ierr=0;
  const unsigned int n_sl = send_list.size();

  IS is;
  VecScatter scatter;

  std::vector<int> idx(n_sl + this->local_size());
  
  for (unsigned int i=0; i<n_sl; i++)
    idx[i] = static_cast<int>(send_list[i]);
  for (unsigned int i = 0; i != this->local_size(); ++i)
    idx[n_sl+i] = i + this->first_local_index();
  
  // Create the index set & scatter object
  if (idx.empty())
    ierr = ISCreateGeneral(libMesh::COMM_WORLD,
                           n_sl+this->local_size(), PETSC_NULL, &is);
  else
    ierr = ISCreateGeneral(libMesh::COMM_WORLD,
                           n_sl+this->local_size(), &idx[0], &is);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);

  ierr = VecScatterCreate(_vec,          is,
                          v_local->_vec, is,
                          &scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  
  // Perform the scatter
#if PETSC_VERSION_LESS_THAN(2,3,3)
         
  ierr = VecScatterBegin(_vec, v_local->_vec, INSERT_VALUES,
                         SCATTER_FORWARD, scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  ierr = VecScatterEnd  (_vec, v_local->_vec, INSERT_VALUES,
                         SCATTER_FORWARD, scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

#else
         
  // API argument order change in PETSc 2.3.3
  ierr = VecScatterBegin(scatter, _vec, v_local->_vec,
                         INSERT_VALUES, SCATTER_FORWARD);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  ierr = VecScatterEnd  (scatter, _vec, v_local->_vec,
                         INSERT_VALUES, SCATTER_FORWARD);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

#endif
         

  // Clean up
  ierr = ISDestroy (is);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  ierr = VecScatterDestroy(scatter);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  // Make sure ghost dofs are up to date
  if (v_local->type() == GHOSTED)
    v_local->close();
}
 

template<typename T > void PetscVector< T >::localize (const unsigned intfirst_local_idx, const unsigned intlast_local_idx, const std::vector< unsigned int > &send_list) [virtual]Updates a local vector with selected values from neighboring processors, as defined by send_list.

Implements NumericVector< T >.

Definition at line 813 of file petsc_vector.C.

References PetscVector< T >::_vec, libMesh::close(), PetscVector< T >::close(), libMesh::COMM_WORLD, MeshTools::Generation::Private::idx(), PetscVector< T >::init(), Utility::iota(), PetscVector< T >::localize(), and libMeshEnums::PARALLEL.

{
  this->_restore_array();

  // Only good for serial vectors.
  // libmesh_assert (this->size() == this->local_size());
  libmesh_assert (last_local_idx > first_local_idx);
  libmesh_assert (send_list.size() <= this->size());
  libmesh_assert (last_local_idx < this->size());
  
  const unsigned int size       = this->size();
  const unsigned int local_size = (last_local_idx - first_local_idx + 1);
  int ierr=0;  
  
  // Don't bother for serial cases
  if ((first_local_idx == 0) &&
      (local_size == size))
    return;
  
  
  // Build a parallel vector, initialize it with the local
  // parts of (*this)
  PetscVector<T> parallel_vec;

  parallel_vec.init (size, local_size, true, PARALLEL);


  // Copy part of *this into the parallel_vec
  {
    IS is;
    VecScatter scatter;

    // Create idx, idx[i] = i+first_local_idx;
    std::vector<int> idx(local_size);
    Utility::iota (idx.begin(), idx.end(), first_local_idx);

    // Create the index set & scatter object
    ierr = ISCreateGeneral(libMesh::COMM_WORLD, local_size, &idx[0], &is); 
           CHKERRABORT(libMesh::COMM_WORLD,ierr);

    ierr = VecScatterCreate(_vec,              is,
                            parallel_vec._vec, is,
                            &scatter);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);

    // Perform the scatter
#if PETSC_VERSION_LESS_THAN(2,3,3)

    ierr = VecScatterBegin(_vec, parallel_vec._vec, INSERT_VALUES,
                           SCATTER_FORWARD, scatter);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
    ierr = VecScatterEnd  (_vec, parallel_vec._vec, INSERT_VALUES,
                           SCATTER_FORWARD, scatter);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);

#else
           
      // API argument order change in PETSc 2.3.3
    ierr = VecScatterBegin(scatter, _vec, parallel_vec._vec,
                           INSERT_VALUES, SCATTER_FORWARD);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
    ierr = VecScatterEnd  (scatter, _vec, parallel_vec._vec,
                           INSERT_VALUES, SCATTER_FORWARD);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);
           
#endif

    // Clean up
    ierr = ISDestroy (is);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
    ierr = VecScatterDestroy(scatter);
           CHKERRABORT(libMesh::COMM_WORLD,ierr);
  }

  // localize like normal
  parallel_vec.close();
  parallel_vec.localize (*this, send_list);
  this->close();
}
 

template<> void PetscVector< Complex >::localize_to_one (std::vector< Complex > &v_local, const unsigned intpid) const

Definition at line 993 of file petsc_vector.C.

References libMesh::COMM_WORLD.

{
  this->_restore_array();

  int ierr=0;
  const int n  = size();
  const int nl = local_size();
  PetscScalar *values;

  
  v_local.resize(n);

  
  for (int i=0; i<n; i++)
    v_local[i] = 0.;
  
  // only one processor
  if (n == nl)
    {      
      ierr = VecGetArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);

      for (int i=0; i<n; i++)
        v_local[i] = static_cast<Complex>(values[i]);

      ierr = VecRestoreArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  // otherwise multiple processors
  else
    {
      unsigned int ioff = this->first_local_index();

      /* in here the local values are stored, acting as send buffer for MPI
       * initialize to zero, since we collect using MPI_SUM
       */
      std::vector<Real> real_local_values(n, 0.);
      std::vector<Real> imag_local_values(n, 0.);

      {
        ierr = VecGetArray (_vec, &values);
               CHKERRABORT(libMesh::COMM_WORLD,ierr);
        
        // provide my local share to the real and imag buffers
        for (int i=0; i<nl; i++)
          {
            real_local_values[i+ioff] = static_cast<Complex>(values[i]).real();
            imag_local_values[i+ioff] = static_cast<Complex>(values[i]).imag();
          }

        ierr = VecRestoreArray (_vec, &values);
               CHKERRABORT(libMesh::COMM_WORLD,ierr);
      }
   
      /* have buffers of the real and imaginary part of v_local.
       * Once MPI_Reduce() collected all the real and imaginary
       * parts in these std::vector<double>, the values can be 
       * copied to v_local
       */
      std::vector<Real> real_v_local(n);
      std::vector<Real> imag_v_local(n);

      // collect entries from other proc's in real_v_local, imag_v_local
      MPI_Reduce (&real_local_values[0], &real_v_local[0], n, 
                  MPI_DOUBLE, MPI_SUM,
                  pid, libMesh::COMM_WORLD);    

      MPI_Reduce (&imag_local_values[0], &imag_v_local[0], n, 
                  MPI_DOUBLE, MPI_SUM,
                  pid, libMesh::COMM_WORLD);    

      // copy real_v_local and imag_v_local to v_local
      for (int i=0; i<n; i++)
        v_local[i] = Complex(real_v_local[i], imag_v_local[i]);
    }  
}
 

template<typename T> void PetscVector< T >::localize_to_one (std::vector< T > &v_local, const unsigned intproc_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 NumericVector< T >.  

template<> void PetscVector< Real >::localize_to_one (std::vector< Real > &v_local, const unsigned intpid) const

Definition at line 936 of file petsc_vector.C.

References libMesh::COMM_WORLD.

{
  this->_restore_array();

  int ierr=0;
  const int n  = size();
  const int nl = local_size();
  PetscScalar *values;

  
  v_local.resize(n);

  
  // only one processor
  if (n == nl)
    {      
      ierr = VecGetArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);

      for (int i=0; i<n; i++)
        v_local[i] = static_cast<Real>(values[i]);

      ierr = VecRestoreArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  // otherwise multiple processors
  else
    {
      unsigned int ioff = this->first_local_index();
      std::vector<Real> local_values (n, 0.);
      
      {
        ierr = VecGetArray (_vec, &values);
               CHKERRABORT(libMesh::COMM_WORLD,ierr);
        
        for (int i=0; i<nl; i++)
          local_values[i+ioff] = static_cast<Real>(values[i]);
        
        ierr = VecRestoreArray (_vec, &values);
               CHKERRABORT(libMesh::COMM_WORLD,ierr);
      }
      

      MPI_Reduce (&local_values[0], &v_local[0], n, MPI_REAL, MPI_SUM,
                  pid, libMesh::COMM_WORLD);
    }
}
 

template<typename T > unsigned int PetscVector< T >::map_global_to_local_index (const unsigned inti) const [inline]Maps the global index i to the corresponding global index. If the index is not a ghost cell, this is done by subtraction the number of the first local index. If it is a ghost cell, it has to be looked up in the map.

Definition at line 1024 of file petsc_vector.h.

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

{
  libmesh_assert (this->initialized());

  int ierr=0, petsc_first=0, petsc_last=0;
  ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
  CHKERRABORT(libMesh::COMM_WORLD,ierr);
  const unsigned int first = static_cast<unsigned int>(petsc_first);
  const unsigned int last = static_cast<unsigned int>(petsc_last);

  if((i>=first) && (i<last))
    {
      return i-first;
    }

  GlobalToLocalMap::const_iterator it = _global_to_local_map.find(i);
  libmesh_assert (it!=_global_to_local_map.end());
  return it->second+last-first;
}
 

template<typename T > Real PetscVector< T >::max () const [inline, virtual]Returns:

the maximum element in the vector. In case of complex numbers, this returns the maximum Real part.

Implements NumericVector< T >.

Definition at line 1110 of file petsc_vector.h.

References libMesh::COMM_WORLD, and std::max().

{
  this->_restore_array();

  int index=0, ierr=0;
  PetscReal max=0.;

  ierr = VecMax (_vec, &index, &max);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  // this return value is correct: VecMax returns a PetscReal
  return static_cast<Real>(max);
}
 

template<typename T > Real PetscVector< T >::min () const [inline, virtual]Returns:

the minimum element in the vector. In case of complex numbers, this returns the minimum Real part.

Implements NumericVector< T >.

Definition at line 1092 of file petsc_vector.h.

References libMesh::COMM_WORLD, and std::min().

{
  this->_restore_array();

  int index=0, ierr=0;
  PetscReal min=0.;

  ierr = VecMin (_vec, &index, &min);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  // this return value is correct: VecMin returns a PetscReal
  return static_cast<Real>(min);
}
 

static unsigned int ReferenceCounter::n_objects () [inline, static, inherited]Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 76 of file reference_counter.h.

References ReferenceCounter::_n_objects.

Referenced by System::read_serialized_blocked_dof_objects(), and System::write_serialized_blocked_dof_objects().

  { return _n_objects; }
 

template<typename T > T PetscVector< T >::operator() (const unsigned inti) const [inline, virtual]Access components, returns U(i).

Implements NumericVector< T >.

Definition at line 1048 of file petsc_vector.h.

References libMeshEnums::GHOSTED.

{
  this->_get_array();

  const unsigned int local_index = this->map_global_to_local_index(i);

#ifndef NDEBUG
  if(this->type() == GHOSTED)
    {
      libmesh_assert(local_index<_local_size);
    }
#endif
  
  return static_cast<T>(_values[local_index]);
}
 

template<typename T> NumericVector<T>& NumericVector< T >::operator*= (const Ta) [inline, inherited]Multiplication operator. Equivalent to U.scale(a)

Definition at line 340 of file numeric_vector.h.

{ this->scale(a); return *this; }
 

template<typename T > NumericVector< T > & PetscVector< T >::operator+= (const NumericVector< T > &V) [virtual]Addition operator. Fast equivalent to U.add(1, V).

Implements NumericVector< T >.

Definition at line 121 of file petsc_vector.C.

References libMesh::closed().

{
  this->_restore_array();
  libmesh_assert(this->closed());
  
  this->add(1., v);
  
  return *this;
}
 

template<typename T > NumericVector< T > & PetscVector< T >::operator-= (const NumericVector< T > &V) [virtual]Subtraction operator. Fast equivalent to U.add(-1, V).

Implements NumericVector< T >.

Definition at line 135 of file petsc_vector.C.

References libMesh::closed().

{
  this->_restore_array();
  libmesh_assert(this->closed());
  
  this->add(-1., v);
  
  return *this;
}
 

template<typename T> NumericVector<T>& NumericVector< T >::operator/= (const Ta) [inline, inherited]Division operator. Equivalent to U.scale(1./a)

Definition at line 346 of file numeric_vector.h.

{ this->scale(1./a); return *this; }
 

template<typename T > NumericVector< T > & PetscVector< T >::operator= (const NumericVector< T > &V) [virtual]$U = V$: copy all components.

Implements NumericVector< T >.

Definition at line 559 of file petsc_vector.C.

{
  // Make sure the NumericVector passed in is really a PetscVector
  const PetscVector<T>* v = libmesh_cast_ptr<const PetscVector<T>*>(&v_in);

  *this = *v;
  
  return *this;
}
 

template<typename T > PetscVector< T > & PetscVector< T >::operator= (const PetscVector< T > &V)$U = V$: copy all components.

Definition at line 573 of file petsc_vector.C.

References PetscVector< T >::_global_to_local_map, PetscVector< T >::_restore_array(), NumericVector< T >::_type, PetscVector< T >::_vec, libMesh::COMM_WORLD, libMeshEnums::GHOSTED, PetscVector< T >::local_size(), and PetscVector< T >::size().

{
  this->_restore_array();
  v._restore_array();
  libmesh_assert (this->_type == v._type);
  libmesh_assert (this->size() == v.size());
  libmesh_assert (this->local_size() == v.local_size());
  libmesh_assert (this->_global_to_local_map ==
                  v._global_to_local_map);

  if (v.size() != 0)
    {
      int ierr = 0;
      if(this->type() != GHOSTED)
        {
          ierr = VecCopy (v._vec, this->_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
        }
      else
        {
          Vec loc_vec;
          Vec v_loc_vec;
          ierr = VecGhostGetLocalForm (_vec,&loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          ierr = VecGhostGetLocalForm (v._vec,&v_loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);

          ierr = VecCopy (v_loc_vec, loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);

          ierr = VecGhostRestoreLocalForm (v._vec,&v_loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
        }
    }
  
  return *this;
}
 

template<typename T > NumericVector< T > & PetscVector< T >::operator= (const std::vector< T > &v) [virtual]$U = V$: copy all components.

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

Case 2: The vector is the same size as our local piece. Insert directly to the local piece.

Implements NumericVector< T >.

Definition at line 617 of file petsc_vector.C.

References libMesh::close(), libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

{
  this->_restore_array();

  const unsigned int nl   = this->local_size();
  const unsigned int ioff = this->first_local_index();
  int ierr=0;
  PetscScalar* values;
      
  if (this->size() == v.size())
    {
      ierr = VecGetArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);

      for (unsigned int i=0; i<nl; i++)
        values[i] =  static_cast<PetscScalar>(v[i+ioff]);
      
      ierr = VecRestoreArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  else
    {
      libmesh_assert (this->local_size() == v.size());

      ierr = VecGetArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);

      for (unsigned int i=0; i<nl; i++)
        values[i] = static_cast<PetscScalar>(v[i]);
      
      ierr = VecRestoreArray (_vec, &values);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  // Make sure ghost dofs are up to date
  if (this->type() == GHOSTED)
    this->close();

  return *this;
}
 

template<typename T > NumericVector< T > & PetscVector< T >::operator= (const Ts) [virtual]Change the dimension to that of the vector V. The same applies as for the other init function.

The elements of V are not copied, i.e. this function is the same as calling init(V.size(),fast). $U(0-N) = s$: fill all components.

Implements NumericVector< T >.

Definition at line 509 of file petsc_vector.C.

References libMesh::closed(), libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

{
  this->_restore_array();
  libmesh_assert(this->closed());

  int ierr = 0;
  PetscScalar s = static_cast<PetscScalar>(s_in);

  if (this->size() != 0)
    {
      if(this->type() != GHOSTED)
        {
#if PETSC_VERSION_LESS_THAN(2,3,0)
          // 2.2.x & earlier style
          ierr = VecSet(&s, _vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
          // 2.3.x & later style         
          ierr = VecSet(_vec, s);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif
        }
      else
        {
          Vec loc_vec;
          ierr = VecGhostGetLocalForm (_vec,&loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);

#if PETSC_VERSION_LESS_THAN(2,3,0)
          // 2.2.x & earlier style
          ierr = VecSet(&s, loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
          // 2.3.x & later style         
          ierr = VecSet(loc_vec, s);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif

          ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
        }
    }
  
  return *this;
}
 

template<typename T > void PetscVector< T >::pointwise_mult (const NumericVector< T > &vec1, const NumericVector< T > &vec2) [virtual]Computes the pointwise (i.e. component-wise) product of vec1 and vec2 and stores the result in *this.

Implements NumericVector< T >.

Definition at line 1077 of file petsc_vector.C.

References libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

{
  this->_restore_array();

  int ierr = 0;

  // Convert arguments to PetscVector*.
  const PetscVector<T>* vec1_petsc = libmesh_cast_ptr<const PetscVector<T>*>(&vec1);
  const PetscVector<T>* vec2_petsc = libmesh_cast_ptr<const PetscVector<T>*>(&vec2);

  // Call PETSc function.

#if PETSC_VERSION_LESS_THAN(2,3,1)

  std::cout << 'This method has been developed with PETSc 2.3.1.  '
            << 'No one has made it backwards compatible with older '
            << 'versions of PETSc so far; however, it might work '
            << 'without any change with some older version.' << std::endl;
  libmesh_error();

#else
  
  if(this->type() != GHOSTED)
    {
      ierr = VecPointwiseMult(this->vec(),
                              const_cast<PetscVector<T>*>(vec1_petsc)->vec(),
                              const_cast<PetscVector<T>*>(vec2_petsc)->vec());
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }
  else
    {
          Vec loc_vec;
          Vec v1_loc_vec;
          Vec v2_loc_vec;
          ierr = VecGhostGetLocalForm (_vec,&loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          ierr = VecGhostGetLocalForm (const_cast<PetscVector<T>*>(vec1_petsc)->vec(),&v1_loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          ierr = VecGhostGetLocalForm (const_cast<PetscVector<T>*>(vec2_petsc)->vec(),&v2_loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);

          ierr = VecPointwiseMult(loc_vec,v1_loc_vec,v2_loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);

          ierr = VecGhostRestoreLocalForm (const_cast<PetscVector<T>*>(vec1_petsc)->vec(),&v1_loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          ierr = VecGhostRestoreLocalForm (const_cast<PetscVector<T>*>(vec2_petsc)->vec(),&v2_loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
          ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
          CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }
      
#endif

}
 

template<> void NumericVector< Complex >::print (std::ostream &os) const [inline, inherited]

Definition at line 731 of file numeric_vector.h.

References libMesh::initialized().

{
  libmesh_assert (this->initialized());
  os << 'Size	global =  ' << this->size()
     << '		local =  ' << this->local_size() << std::endl;
  
  // std::complex<>::operator<<() is defined, but use this form
  os << '#	Real part		Imaginary part' << std::endl;
  for (unsigned int i=this->first_local_index(); i<this->last_local_index(); i++)
    os << i << '	' 
       << (*this)(i).real() << '		' 
       << (*this)(i).imag() << std::endl;
}
 

template<typename T > void NumericVector< T >::print (std::ostream &os = std::cout) const [inline, virtual, inherited]Prints the local contents of the vector to the screen.

Definition at line 749 of file numeric_vector.h.

References libMesh::initialized().

{
  libmesh_assert (this->initialized());
  os << 'Size	global =  ' << this->size()
     << '		local =  ' << this->local_size() << std::endl;

  os << '#	Value' << std::endl;
  for (unsigned int i=this->first_local_index(); i<this->last_local_index(); i++)
    os << i << '	' << (*this)(i) << std::endl;
}
 

template<typename T > void NumericVector< T >::print_global (std::ostream &os = std::cout) const [inline, virtual, inherited]Prints the global contents of the vector to the screen.

Definition at line 786 of file numeric_vector.h.

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

{
  libmesh_assert (this->initialized());

  std::vector<T> v(this->size());
  this->localize(v);

  // Right now we only want one copy of the output
  if (libMesh::processor_id())
    return;

  os << 'Size	global =  ' << this->size() << std::endl;
  os << '#	Value' << std::endl;
  for (unsigned int i=0; i!=v.size(); i++)
    os << i << '	' << v[i] << std::endl;
}
 

template<> void NumericVector< Complex >::print_global (std::ostream &os) const [inline, inherited]

Definition at line 764 of file numeric_vector.h.

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

{
  libmesh_assert (this->initialized());

  std::vector<Complex> v(this->size());
  this->localize(v);

  // Right now we only want one copy of the output
  if (libMesh::processor_id())
    return;
  
  os << 'Size	global =  ' << this->size() << std::endl;
  os << '#	Real part		Imaginary part' << std::endl;
  for (unsigned int i=0; i!=v.size(); i++)
    os << i << '	' 
       << v[i].real() << '		' 
       << v[i].imag() << std::endl;
}
 

void ReferenceCounter::print_info () [static, inherited]Prints the reference information to std::cout.

Definition at line 83 of file reference_counter.C.

References ReferenceCounter::get_info().

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
  
  std::cout << ReferenceCounter::get_info();
  
#endif
}
 

template<typename T > void PetscVector< T >::print_matlab (const std::stringname = 'NULL') const [virtual]Print the contents of the vector in Matlab format. Optionally prints the matrix to the file named name. If name is not specified it is dumped to the screen.

Create an ASCII file containing the matrix if a filename was provided.

Otherwise the matrix will be dumped to the screen.

Destroy the viewer.

Reimplemented from NumericVector< T >.

Definition at line 1137 of file petsc_vector.C.

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

{
  this->_restore_array();
  libmesh_assert (this->closed());
  
  int ierr=0; 
  PetscViewer petsc_viewer;


  ierr = PetscViewerCreate (libMesh::COMM_WORLD,
                            &petsc_viewer);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  if (name != 'NULL')
    {
      ierr = PetscViewerASCIIOpen( libMesh::COMM_WORLD,
                                   name.c_str(),
                                   &petsc_viewer);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
      
      ierr = PetscViewerSetFormat (petsc_viewer,
                                   PETSC_VIEWER_ASCII_MATLAB);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
      ierr = VecView (_vec, petsc_viewer);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }

  else
    {
      ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
                                   PETSC_VIEWER_ASCII_MATLAB);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
      ierr = VecView (_vec, PETSC_VIEWER_STDOUT_WORLD);
             CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }


  ierr = PetscViewerDestroy (petsc_viewer);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
}
 

template<typename T > void PetscVector< T >::scale (const Tfactor) [virtual]Scale each element of the vector by the given factor.

Implements NumericVector< T >.

Definition at line 417 of file petsc_vector.C.

References libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

{
  this->_restore_array();

  int ierr = 0;
  PetscScalar factor = static_cast<PetscScalar>(factor_in);
  
  if(this->type() != GHOSTED)
    {
#if PETSC_VERSION_LESS_THAN(2,3,0)
      // 2.2.x & earlier style
      ierr = VecScale(&factor, _vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
      // 2.3.x & later style     
      ierr = VecScale(_vec, factor);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif
    }
  else
    {
      Vec loc_vec;
      ierr = VecGhostGetLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);

#if PETSC_VERSION_LESS_THAN(2,3,0)
      // 2.2.x & earlier style
      ierr = VecScale(&factor, loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
      // 2.3.x & later style     
      ierr = VecScale(loc_vec, factor);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif

      ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }
}
 

template<typename T > void PetscVector< T >::set (const unsigned inti, const Tvalue) [virtual]v(i) = value

Implements NumericVector< T >.

Definition at line 148 of file petsc_vector.C.

References libMesh::COMM_WORLD.

{
  this->_restore_array();
  libmesh_assert(i<size());
  
  int ierr=0;
  int i_val = static_cast<int>(i);
  PetscScalar petsc_value = static_cast<PetscScalar>(value);

  ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, INSERT_VALUES);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  this->_is_closed = false;
}
 

template<typename T > unsigned int PetscVector< T >::size () const [inline, virtual]Returns:

dimension of the vector. This function was formerly called n(), but was renamed to get the PetscVector<T> class closer to the C++ standard library's std::vector container.

Implements NumericVector< T >.

Definition at line 957 of file petsc_vector.h.

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

Referenced by PetscVector< T >::add(), PetscVector< T >::init(), PetscVector< T >::localize(), and PetscVector< T >::operator=().

{
  libmesh_assert (this->initialized());
  
  int ierr=0, petsc_size=0;
  
  if (!this->initialized())
    return 0;
  
  ierr = VecGetSize(_vec, &petsc_size);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);

  return static_cast<unsigned int>(petsc_size);
}
 

template<class T > Real NumericVector< T >::subset_l1_norm (const std::set< unsigned int > &indices) [virtual, inherited]Returns:

the $l_1$-norm of the vector, i.e. the sum of the absolute values for the specified entries in the vector.

Note that the indices must necessary live on this processor.

Definition at line 193 of file numeric_vector.C.

Referenced by System::discrete_var_norm().

{
  NumericVector<T> & v = *this;
  
  std::set<unsigned int>::iterator it = indices.begin();
  const std::set<unsigned int>::iterator it_end = indices.end();

  Real norm = 0;
  
  for(; it!=it_end; ++it)
    norm += std::abs(v(*it));

  Parallel::sum(norm);

  return norm;
}
 

template<class T > Real NumericVector< T >::subset_l2_norm (const std::set< unsigned int > &indices) [virtual, inherited]Returns:

the $l_2$-norm of the vector, i.e. the square root of the sum of the squares of the elements for the specified entries in the vector.

Note that the indices must necessary live on this processor.

Definition at line 211 of file numeric_vector.C.

References libmesh_norm().

Referenced by System::discrete_var_norm().

{
  NumericVector<T> & v = *this;

  std::set<unsigned int>::iterator it = indices.begin();
  const std::set<unsigned int>::iterator it_end = indices.end();

  Real norm = 0;
  
  for(; it!=it_end; ++it)
    norm += libmesh_norm(v(*it));

  Parallel::sum(norm);

  return std::sqrt(norm);
}
 

template<class T > Real NumericVector< T >::subset_linfty_norm (const std::set< unsigned int > &indices) [virtual, inherited]Returns:

the maximum absolute value of the specified entries of this vector, which is the $l_infty$-norm of a vector.

Note that the indices must necessary live on this processor.

Definition at line 229 of file numeric_vector.C.

References std::max().

Referenced by System::discrete_var_norm().

{
  NumericVector<T> & v = *this;

  std::set<unsigned int>::iterator it = indices.begin();
  const std::set<unsigned int>::iterator it_end = indices.end();

  Real norm = 0;
  
  for(; it!=it_end; ++it)
    {
      Real value = std::abs(v(*it));
      if(value > norm)
        norm = value;
    }

  Parallel::max(norm);

  return norm;
}
 

template<typename T > T PetscVector< T >::sum () const [virtual]Returns:

the sum of values in a vector

Implements NumericVector< T >.

Definition at line 51 of file petsc_vector.C.

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

{
  this->_restore_array();
  libmesh_assert(this->closed());
  
  int ierr=0;
  PetscScalar value=0.;
  
  ierr = VecSum (_vec, &value);
         CHKERRABORT(libMesh::COMM_WORLD,ierr);
  
  return static_cast<T>(value);
}
 

template<typename T > void PetscVector< T >::swap (NumericVector< T > &v) [inline, virtual]Swaps the raw PETSc vector context pointers.

Reimplemented from NumericVector< T >.

Definition at line 1128 of file petsc_vector.h.

References PetscVector< T >::_array_is_present, PetscVector< T >::_destroy_vec_on_exit, PetscVector< T >::_global_to_local_map, PetscVector< T >::_local_form, PetscVector< T >::_values, PetscVector< T >::_vec, and NumericVector< T >::swap().

{
  NumericVector<T>::swap(other);

  PetscVector<T>& v = libmesh_cast_ref<PetscVector<T>&>(other);

  std::swap(_vec, v._vec);
  std::swap(_destroy_vec_on_exit, v._destroy_vec_on_exit);
  std::swap(_global_to_local_map, v._global_to_local_map);
  std::swap(_array_is_present, v._array_is_present);
  std::swap(_local_form, v._local_form);
  std::swap(_values, v._values);
}
 

template<typename T> ParallelType& NumericVector< T >::type () [inline, inherited]Returns:

the type (SERIAL, PARALLEL, GHOSTED) of the vector.

Definition at line 119 of file numeric_vector.h.

{ return _type; }
 

template<typename T> ParallelType NumericVector< T >::type () const [inline, inherited]Returns:

the type (SERIAL, PARALLEL, GHOSTED) of the vector.

Definition at line 114 of file numeric_vector.h.

Referenced by DistributedVector< T >::DistributedVector(), DofMap::enforce_constraints_exactly(), EpetraVector< T >::EpetraVector(), EpetraVector< T >::init(), LaspackVector< T >::init(), DistributedVector< T >::init(), PetscVector< T >::localize(), PetscVector< T >::PetscVector(), System::project_vector(), and System::read_serialized_vector().

{ return _type; }
 

template<typename T> Vec PetscVector< T >::vec () [inline]Returns the raw PETSc vector context pointer. Note this is generally not required in user-level code. Just don't do anything crazy like calling VecDestroy()!

Definition at line 486 of file petsc_vector.h.

References PetscVector< T >::_vec.

Referenced by PetscPreconditioner< T >::apply(), PetscMatrix< T >::get_diagonal(), SlepcEigenSolver< T >::get_eigenpair(), PetscLinearSolver< T >::solve(), and PetscDiffSolver::solve().

{ libmesh_assert (_vec != NULL); return _vec; }
 

template<typename T > void PetscVector< T >::zero () [inline, virtual]Set all entries to zero. Equivalent to v = 0, but more obvious and faster.

Implements NumericVector< T >.

Definition at line 897 of file petsc_vector.h.

References libMesh::COMM_WORLD, and libMeshEnums::GHOSTED.

{
  this->_restore_array();
  
  int ierr=0;

  PetscScalar z=0.;

  if(this->type() != GHOSTED)
    {
#if PETSC_VERSION_LESS_THAN(2,3,0)  
      // 2.2.x & earlier style
      ierr = VecSet (&z, _vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
      // 2.3.x & newer
      ierr = VecSet (_vec, z);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif
    }
  else
    {
      /* Vectors that include ghost values require a special
         handling.  */
      Vec loc_vec;
      ierr = VecGhostGetLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#if PETSC_VERSION_LESS_THAN(2,3,0)  
      // 2.2.x & earlier style
      ierr = VecSet (&z, loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#else
      // 2.3.x & newer
      ierr = VecSet (loc_vec, z);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
#endif
      ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
      CHKERRABORT(libMesh::COMM_WORLD,ierr);
    }
}
 

Friends And Related Function Documentation

 

template<typename T> std::ostream& operator<< (std::ostream &os, const NumericVector< T > &v) [friend, inherited]Same as above but allows you to use stream syntax.

Definition at line 537 of file numeric_vector.h.

  {
    v.print_global(os);
    return os;
  }
 

Member Data Documentation

 

template<typename T> bool PetscVector< T >::_array_is_present [mutable, private]If true, the actual Petsc array of the values of the vector is currently accessible. That means that the members _local_form and _values are valid.

Definition at line 503 of file petsc_vector.h.

Referenced by PetscVector< T >::swap().  

ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited]Actually holds the data.

Definition at line 110 of file reference_counter.h.

Referenced by ReferenceCounter::get_info(), ReferenceCounter::increment_constructor_count(), and ReferenceCounter::increment_destructor_count().  

template<typename T> bool PetscVector< T >::_destroy_vec_on_exit [private]This boolean value should only be set to false for the constructor which takes a PETSc Vec object.

Definition at line 554 of file petsc_vector.h.

Referenced by PetscVector< T >::swap().  

template<typename T> GlobalToLocalMap PetscVector< T >::_global_to_local_map [private]Map that maps global to local ghost cells (will be empty if not in ghost cell mode).

Definition at line 548 of file petsc_vector.h.

Referenced by PetscVector< T >::init(), PetscVector< T >::operator=(), PetscVector< T >::PetscVector(), and PetscVector< T >::swap().  

template<typename T> bool NumericVector< T >::_is_closed [protected, inherited]Flag to see if the Numeric assemble routines have been called yet

Definition at line 582 of file numeric_vector.h.

Referenced by NumericVector< Number >::closed(), PetscVector< T >::init(), DistributedVector< T >::localize(), LaspackVector< T >::operator=(), DistributedVector< T >::operator=(), PetscVector< T >::PetscVector(), and NumericVector< T >::swap().  

template<typename T> bool NumericVector< T >::_is_initialized [protected, inherited]Flag to tell if init has been called yet

Definition at line 588 of file numeric_vector.h.

Referenced by PetscVector< T >::create_subvector(), EpetraVector< T >::EpetraVector(), PetscVector< T >::init(), NumericVector< Number >::initialized(), DistributedVector< T >::localize(), DistributedVector< T >::operator=(), PetscVector< T >::PetscVector(), and NumericVector< T >::swap().  

template<typename T> Vec PetscVector< T >::_local_form [mutable, private]Petsc vector datatype to hold the local form of a ghosted vector. The contents of this field are only valid if the vector is ghosted and _array_is_present is true.

Definition at line 519 of file petsc_vector.h.

Referenced by PetscVector< T >::swap().  

template<typename T> unsigned int PetscVector< T >::_local_size [mutable, private]Size of the local form, for being used in assertations. The contents of this field are only valid if the vector is ghosted and _array_is_present is true.

Definition at line 511 of file petsc_vector.h.  

Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited]Mutual exclusion object to enable thread-safe reference counting.

Definition at line 123 of file reference_counter.h.  

Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited]The number of objects. Print the reference count information when the number returns to 0.

Definition at line 118 of file reference_counter.h.

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

template<typename T> ParallelType NumericVector< T >::_type [protected, inherited]Type of vector

Definition at line 593 of file numeric_vector.h.

Referenced by DistributedVector< T >::DistributedVector(), EpetraVector< T >::EpetraVector(), PetscVector< T >::init(), PetscVector< T >::operator=(), PetscVector< T >::PetscVector(), NumericVector< T >::swap(), and NumericVector< Number >::type().  

template<typename T> PetscScalar* PetscVector< T >::_values [mutable, private]Pointer to the actual Petsc array of the values of the vector. This pointer is only valid if _array_is_present is true.

Definition at line 525 of file petsc_vector.h.

Referenced by PetscVector< T >::swap().  

template<typename T> Vec PetscVector< T >::_vec [private]Actual Petsc vector datatype to hold vector entries

Definition at line 496 of file petsc_vector.h.

Referenced by PetscVector< T >::add(), PetscVector< T >::add_vector(), PetscVector< T >::create_subvector(), PetscVector< T >::dot(), PetscVector< T >::init(), PetscVector< T >::localize(), PetscVector< T >::operator=(), PetscVector< T >::PetscVector(), PetscVector< T >::swap(), and PetscVector< T >::vec().

 

Author

Generated automatically by Doxygen for libMesh from the source code.


 

Index

NAME
SYNOPSIS
Public Member Functions
Static Public Member Functions
Protected Types
Protected Member Functions
Protected Attributes
Static Protected Attributes
Private Types
Private Member Functions
Private Attributes
Friends
Detailed Description
template<typename T> class PetscVector< T >
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited]Data structure to log the information. The log is identified by the class name.
template<typename T> typedef std::map<unsigned int,unsigned int> PetscVector< T >::GlobalToLocalMap [private]Type for map that maps global to local ghost cells.
Constructor & Destructor Documentation
template<typename T > PetscVector< T >::PetscVector (const ParallelTypetype = AUTOMATIC) [inline, explicit]Dummy-Constructor. Dimension=0
template<typename T > PetscVector< T >::PetscVector (const unsigned intn, const ParallelTypetype = AUTOMATIC) [inline, explicit]Constructor. Set dimension to n and initialize all elements with zero.
template<typename T > PetscVector< T >::PetscVector (const unsigned intn, const unsigned intn_local, const ParallelTypetype = AUTOMATIC) [inline]Constructor. Set local dimension to n_local, the global dimension to n, and initialize all elements with zero.
template<typename T > PetscVector< T >::PetscVector (const unsigned intN, const unsigned intn_local, const std::vector< unsigned int > &ghost, const ParallelTypetype = 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.
template<typename T > PetscVector< T >::PetscVector (Vecv) [inline]Constructor. Creates a PetscVector assuming you already have a valid PETSc Vec object. In this case, v is NOT destroyed by the PetscVector constructor when this object goes out of scope. This allows ownership of v to remain with the original creator, and to simply provide additional functionality with the PetscVector.
template<typename T > PetscVector< T >::~PetscVector () [inline]Destructor, deallocates memory. Made virtual to allow for derived classes to behave properly.
Member Function Documentation
template<typename T > void PetscVector< T >::_get_array (void) const [inline, private]Queries the array (and the local form if the vector is ghosted) from Petsc.
template<typename T > void PetscVector< T >::_restore_array (void) const [inline, private]Restores the array (and the local form if the vector is ghosted) to Petsc.
template<typename T > void PetscVector< T >::abs () [virtual]v = abs(v)... that is, each entry in v is replaced by its absolute value.
template<typename T > void PetscVector< T >::add (const unsigned inti, const Tvalue) [virtual]v(i) += value
template<typename T > void PetscVector< T >::add (const Ta, const NumericVector< T > &v) [virtual]$ U+=a*V $ . Simple vector addition, equal to the operator +=.
template<typename T > void PetscVector< T >::add (const Ts) [virtual]$U(0-LIBMESH_DIM)+=s$. Addition of s to all components. Note that s is a scalar and not a vector.
template<typename T > void PetscVector< T >::add (const NumericVector< T > &V) [virtual]$ U+=V $ . Simple vector addition, equal to the operator +=.
template<typename T > void PetscVector< T >::add_vector (const std::vector< T > &v, const std::vector< unsigned int > &dof_indices) [virtual]$ U+=v $ where v is a std::vector<T> and you want to specify WHERE to add it
template<typename T > void PetscVector< T >::add_vector (const NumericVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$ U+=V $ where U and V are type NumericVector<T> and you want to specify WHERE to add the NumericVector<T> V
template<typename T > void PetscVector< T >::add_vector (const NumericVector< T > &V, const SparseMatrix< T > &A) [virtual]$U+=A*V$, add the product of a SparseMatrix A and a NumericVector V to this, where this=U.
template<typename T> void NumericVector< T >::add_vector (const NumericVector< T > &v, const ShellMatrix< T > &a) [inherited]$U+=A*V$, add the product of a ShellMatrix A and a NumericVector V to this, where this=U.
template<typename T > void PetscVector< T >::add_vector (const DenseVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$U+=V $ where U and V are type DenseVector<T> and you want to specify WHERE to add the DenseVector<T> V
template<typename T > AutoPtr< NumericVector< T > > NumericVector< T >::build (const SolverPackagesolver_package = libMesh::default_solver_package()) [static, inherited]Builds a NumericVector using the linear solver package specified by solver_package
template<typename T > void PetscVector< T >::clear () [inline, virtual]Returns:
template<typename T > AutoPtr< NumericVector< T > > PetscVector< T >::clone () const [inline, virtual]Creates a copy of this vector and returns it in an AutoPtr.
template<typename T > void PetscVector< T >::close () [inline, virtual]Call the assemble functions
template<typename T> virtual bool NumericVector< T >::closed () const [inline, virtual, inherited]Returns:
template<typename T> virtual int NumericVector< T >::compare (const NumericVector< T > &other_vector, const Realthreshold = TOLERANCE) const [virtual, inherited]Returns:
template<> int NumericVector< float >::compare (const NumericVector< float > &other_vector, const Realthreshold) const [inherited]
template<> int NumericVector< double >::compare (const NumericVector< double > &other_vector, const Realthreshold) const [inherited]
template<> int NumericVector< long double >::compare (const NumericVector< long double > &other_vector, const Realthreshold) const [inherited]
template<> int NumericVector< Complex >::compare (const NumericVector< Complex > &other_vector, const Realthreshold) const [inherited]
template<typename T > void PetscVector< T >::create_subvector (NumericVector< T > &subvector, const std::vector< unsigned int > &rows) const [virtual]Creates a 'subvector' from this vector using the rows indices of the 'rows' array.
template<typename T > T PetscVector< T >::dot (const NumericVector< T > &V) const [virtual]Computes the dot product, p = U.V
template<typename T> virtual T NumericVector< T >::el (const unsigned inti) const [inline, virtual, inherited]Returns:
template<typename T > unsigned int PetscVector< T >::first_local_index () const [inline, virtual]Returns:
template<typename T > void PetscVector< T >::get (const std::vector< unsigned int > &index, std::vector< T > &values) const [inline, virtual]Access multiple components at once. Overloaded method that should be faster (probably much faster) than calling operator() individually for each index.
std::string ReferenceCounter::get_info () [static, inherited]Gets a string containing the reference information.
void ReferenceCounter::increment_constructor_count (const std::string &name) [inline, protected, inherited]Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
void ReferenceCounter::increment_destructor_count (const std::string &name) [inline, protected, inherited]Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
template<typename T > void PetscVector< T >::init (const NumericVector< T > &other, const boolfast = false) [inline, virtual]Creates a vector that has the same dimension and storage type as other, including ghost dofs.
template<typename T > void PetscVector< T >::init (const unsigned intN, const unsigned intn_local, const boolfast = false, const ParallelTypetype = AUTOMATIC) [inline, virtual]Change the dimension of the vector to N. The reserved memory for this vector remains unchanged if possible, to make things faster, but this may waste some memory, so take this in the back of your head. However, if N==0 all memory is freed, i.e. if you want to resize the vector and release the memory not needed, you have to first call init(0) and then init(N). This cited behaviour is analogous to that of the STL containers.
template<typename T > void PetscVector< T >::init (const unsigned intN, const boolfast = false, const ParallelTypetype = AUTOMATIC) [inline, virtual]call init with n_local = N,
template<typename T > void PetscVector< T >::init (const unsigned intn, const unsigned intn_local, const std::vector< unsigned int > &ghost, const boolfast = false, const ParallelTypetype = AUTOMATIC) [inline, virtual]Create a vector that holds tha local indices plus those specified in the ghost argument.
template<typename T> virtual bool NumericVector< T >::initialized () const [inline, virtual, inherited]Returns:
template<typename T > void PetscVector< T >::insert (const std::vector< T > &v, const std::vector< unsigned int > &dof_indices) [virtual]$ U=v $ where v is a std::vector<T> and you want to specify WHERE to insert it
template<typename T > void PetscVector< T >::insert (const NumericVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$U=V$, where U and V are type NumericVector<T> and you want to specify WHERE to insert the NumericVector<T> V
template<typename T > void PetscVector< T >::insert (const DenseSubVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$ U=V $ where V is type DenseSubVector<T> and you want to specify WHERE to insert it
template<typename T > void PetscVector< T >::insert (const DenseVector< T > &V, const std::vector< unsigned int > &dof_indices) [virtual]$ U=V $ where V is type DenseVector<T> and you want to specify WHERE to insert it
template<typename T > Real PetscVector< T >::l1_norm () const [virtual]Returns:
template<typename T > Real PetscVector< T >::l2_norm () const [virtual]Returns:
template<typename T > unsigned int PetscVector< T >::last_local_index () const [inline, virtual]Returns:
template<typename T > Real PetscVector< T >::linfty_norm () const [virtual]Returns:
template<typename T > unsigned int PetscVector< T >::local_size () const [inline, virtual]Returns:
template<typename T > void PetscVector< T >::localize (std::vector< T > &v_local) const [virtual]Creates a copy of the global vector in the local vector v_local.
template<typename T > void PetscVector< T >::localize (NumericVector< T > &v_local) const [virtual]Same, but fills a NumericVector<T> instead of a std::vector.
template<typename T > void PetscVector< T >::localize (NumericVector< T > &v_local, const std::vector< unsigned int > &send_list) const [virtual]Creates a local vector v_local containing only information relevant to this processor, as defined by the send_list.
template<typename T > void PetscVector< T >::localize (const unsigned intfirst_local_idx, const unsigned intlast_local_idx, const std::vector< unsigned int > &send_list) [virtual]Updates a local vector with selected values from neighboring processors, as defined by send_list.
template<> void PetscVector< Complex >::localize_to_one (std::vector< Complex > &v_local, const unsigned intpid) const
template<typename T> void PetscVector< T >::localize_to_one (std::vector< T > &v_local, const unsigned intproc_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.
template<> void PetscVector< Real >::localize_to_one (std::vector< Real > &v_local, const unsigned intpid) const
template<typename T > unsigned int PetscVector< T >::map_global_to_local_index (const unsigned inti) const [inline]Maps the global index i to the corresponding global index. If the index is not a ghost cell, this is done by subtraction the number of the first local index. If it is a ghost cell, it has to be looked up in the map.
template<typename T > Real PetscVector< T >::max () const [inline, virtual]Returns:
template<typename T > Real PetscVector< T >::min () const [inline, virtual]Returns:
static unsigned int ReferenceCounter::n_objects () [inline, static, inherited]Prints the number of outstanding (created, but not yet destroyed) objects.
template<typename T > T PetscVector< T >::operator() (const unsigned inti) const [inline, virtual]Access components, returns U(i).
template<typename T> NumericVector<T>& NumericVector< T >::operator*= (const Ta) [inline, inherited]Multiplication operator. Equivalent to U.scale(a)
template<typename T > NumericVector< T > & PetscVector< T >::operator+= (const NumericVector< T > &V) [virtual]Addition operator. Fast equivalent to U.add(1, V).
template<typename T > NumericVector< T > & PetscVector< T >::operator-= (const NumericVector< T > &V) [virtual]Subtraction operator. Fast equivalent to U.add(-1, V).
template<typename T> NumericVector<T>& NumericVector< T >::operator/= (const Ta) [inline, inherited]Division operator. Equivalent to U.scale(1./a)
template<typename T > NumericVector< T > & PetscVector< T >::operator= (const NumericVector< T > &V) [virtual]$U = V$: copy all components.
template<typename T > PetscVector< T > & PetscVector< T >::operator= (const PetscVector< T > &V)$U = V$: copy all components.
template<typename T > NumericVector< T > & PetscVector< T >::operator= (const std::vector< T > &v) [virtual]$U = V$: copy all components.
template<typename T > NumericVector< T > & PetscVector< T >::operator= (const Ts) [virtual]Change the dimension to that of the vector V. The same applies as for the other init function.
template<typename T > void PetscVector< T >::pointwise_mult (const NumericVector< T > &vec1, const NumericVector< T > &vec2) [virtual]Computes the pointwise (i.e. component-wise) product of vec1 and vec2 and stores the result in *this.
template<> void NumericVector< Complex >::print (std::ostream &os) const [inline, inherited]
template<typename T > void NumericVector< T >::print (std::ostream &os = std::cout) const [inline, virtual, inherited]Prints the local contents of the vector to the screen.
template<typename T > void NumericVector< T >::print_global (std::ostream &os = std::cout) const [inline, virtual, inherited]Prints the global contents of the vector to the screen.
template<> void NumericVector< Complex >::print_global (std::ostream &os) const [inline, inherited]
void ReferenceCounter::print_info () [static, inherited]Prints the reference information to std::cout.
template<typename T > void PetscVector< T >::print_matlab (const std::stringname = 'NULL') const [virtual]Print the contents of the vector in Matlab format. Optionally prints the matrix to the file named name. If name is not specified it is dumped to the screen.
template<typename T > void PetscVector< T >::scale (const Tfactor) [virtual]Scale each element of the vector by the given factor.
template<typename T > void PetscVector< T >::set (const unsigned inti, const Tvalue) [virtual]v(i) = value
template<typename T > unsigned int PetscVector< T >::size () const [inline, virtual]Returns:
template<class T > Real NumericVector< T >::subset_l1_norm (const std::set< unsigned int > &indices) [virtual, inherited]Returns:
template<class T > Real NumericVector< T >::subset_l2_norm (const std::set< unsigned int > &indices) [virtual, inherited]Returns:
template<class T > Real NumericVector< T >::subset_linfty_norm (const std::set< unsigned int > &indices) [virtual, inherited]Returns:
template<typename T > T PetscVector< T >::sum () const [virtual]Returns:
template<typename T > void PetscVector< T >::swap (NumericVector< T > &v) [inline, virtual]Swaps the raw PETSc vector context pointers.
template<typename T> ParallelType& NumericVector< T >::type () [inline, inherited]Returns:
template<typename T> ParallelType NumericVector< T >::type () const [inline, inherited]Returns:
template<typename T> Vec PetscVector< T >::vec () [inline]Returns the raw PETSc vector context pointer. Note this is generally not required in user-level code. Just don't do anything crazy like calling VecDestroy()!
template<typename T > void PetscVector< T >::zero () [inline, virtual]Set all entries to zero. Equivalent to v = 0, but more obvious and faster.
Friends And Related Function Documentation
template<typename T> std::ostream& operator<< (std::ostream &os, const NumericVector< T > &v) [friend, inherited]Same as above but allows you to use stream syntax.
Member Data Documentation
template<typename T> bool PetscVector< T >::_array_is_present [mutable, private]If true, the actual Petsc array of the values of the vector is currently accessible. That means that the members _local_form and _values are valid.
ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited]Actually holds the data.
template<typename T> bool PetscVector< T >::_destroy_vec_on_exit [private]This boolean value should only be set to false for the constructor which takes a PETSc Vec object.
template<typename T> GlobalToLocalMap PetscVector< T >::_global_to_local_map [private]Map that maps global to local ghost cells (will be empty if not in ghost cell mode).
template<typename T> bool NumericVector< T >::_is_closed [protected, inherited]Flag to see if the Numeric assemble routines have been called yet
template<typename T> bool NumericVector< T >::_is_initialized [protected, inherited]Flag to tell if init has been called yet
template<typename T> Vec PetscVector< T >::_local_form [mutable, private]Petsc vector datatype to hold the local form of a ghosted vector. The contents of this field are only valid if the vector is ghosted and _array_is_present is true.
template<typename T> unsigned int PetscVector< T >::_local_size [mutable, private]Size of the local form, for being used in assertations. The contents of this field are only valid if the vector is ghosted and _array_is_present is true.
Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited]Mutual exclusion object to enable thread-safe reference counting.
Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited]The number of objects. Print the reference count information when the number returns to 0.
template<typename T> ParallelType NumericVector< T >::_type [protected, inherited]Type of vector
template<typename T> PetscScalar* PetscVector< T >::_values [mutable, private]Pointer to the actual Petsc array of the values of the vector. This pointer is only valid if _array_is_present is true.
template<typename T> Vec PetscVector< T >::_vec [private]Actual Petsc vector datatype to hold vector entries
Author

This document was created by man2html, using the manual pages.
Time: 21:52:24 GMT, April 16, 2011