EpetraVector< T > Class Template Reference

#include <trilinos_epetra_vector.h>

Inheritance diagram for EpetraVector< T >:

List of all members.

Public Member Functions

 EpetraVector (const ParallelType type=AUTOMATIC)
 EpetraVector (const unsigned int n, const ParallelType type=AUTOMATIC)
 EpetraVector (const unsigned int n, const unsigned int n_local, const ParallelType type=AUTOMATIC)
 EpetraVector (const unsigned int N, const unsigned int n_local, const std::vector< unsigned int > &ghost, const ParallelType type=AUTOMATIC)
 EpetraVector (Epetra_Vector &v)
 ~EpetraVector ()
void close ()
void clear ()
void zero ()
virtual AutoPtr< NumericVector
< T > > 
zero_clone () const
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)
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)
EpetraVector< T > & operator= (const EpetraVector< T > &V)
NumericVector< T > & operator= (const std::vector< T > &v)
Real min () const
Real max () const
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
operator() (const unsigned int i) 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)
Epetra_Vector * vec ()
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
virtual void get (const std::vector< unsigned int > &index, std::vector< T > &values) 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 Member Functions

int SumIntoGlobalValues (int numIDs, const int *GIDs, const double *values)
int SumIntoGlobalValues (const Epetra_IntSerialDenseVector &GIDs, const Epetra_SerialDenseVector &values)
int ReplaceGlobalValues (int numIDs, const int *GIDs, const double *values)
int ReplaceGlobalValues (const Epetra_IntSerialDenseVector &GIDs, const Epetra_SerialDenseVector &values)
int SumIntoGlobalValues (int numIDs, const int *GIDs, const int *numValuesPerID, const double *values)
int ReplaceGlobalValues (int numIDs, const int *GIDs, const int *numValuesPerID, const double *values)
int GlobalAssemble (Epetra_CombineMode mode=Add)
void setIgnoreNonLocalEntries (bool flag)
void FEoperatorequals (const EpetraVector &source)
int inputValues (int numIDs, const int *GIDs, const double *values, bool accumulate)
int inputValues (int numIDs, const int *GIDs, const int *numValuesPerID, const double *values, bool accumulate)
int inputNonlocalValue (int GID, double value, bool accumulate)
int inputNonlocalValues (int GID, int numValues, const double *values, bool accumulate)
void destroyNonlocalData ()

Private Attributes

Epetra_Vector * _vec
Epetra_Map * _map
bool _destroy_vec_on_exit
int myFirstID_
int myNumIDs_
double * myCoefs_
int * nonlocalIDs_
int * nonlocalElementSize_
int numNonlocalIDs_
int allocatedNonlocalLength_
double ** nonlocalCoefs_
bool ignoreNonLocalEntries_

Friends

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


Detailed Description

template<typename T>
class EpetraVector< T >

Epetra vector. Provides a nice interface to the Trilinos Epetra data structures for parallel vectors.

Author:
Derek R. Gaston, 2008

Definition at line 55 of file trilinos_epetra_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.


Constructor & Destructor Documentation

template<typename T >
EpetraVector< T >::EpetraVector ( const ParallelType  type = AUTOMATIC  )  [inline, explicit]

Dummy-Constructor. Dimension=0

Definition at line 599 of file trilinos_epetra_vector.h.

References NumericVector< T >::_type.

00600 : _destroy_vec_on_exit(true),
00601   myFirstID_(0),
00602   myNumIDs_(0),
00603   myCoefs_(NULL),
00604   nonlocalIDs_(NULL),
00605   nonlocalElementSize_(NULL),
00606   numNonlocalIDs_(0),
00607   allocatedNonlocalLength_(0),
00608   nonlocalCoefs_(NULL),
00609   ignoreNonLocalEntries_(false)
00610 {
00611   this->_type = type;
00612 }

template<typename T >
EpetraVector< T >::EpetraVector ( const unsigned int  n,
const ParallelType  type = AUTOMATIC 
) [inline, explicit]

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

Definition at line 618 of file trilinos_epetra_vector.h.

References EpetraVector< T >::init().

00620 : _destroy_vec_on_exit(true),
00621   myFirstID_(0),
00622   myNumIDs_(0),
00623   myCoefs_(NULL),
00624   nonlocalIDs_(NULL),
00625   nonlocalElementSize_(NULL),
00626   numNonlocalIDs_(0),
00627   allocatedNonlocalLength_(0),
00628   nonlocalCoefs_(NULL),
00629   ignoreNonLocalEntries_(false)
00630 
00631 {
00632   this->init(n, n, false, type);
00633 }

template<typename T >
EpetraVector< T >::EpetraVector ( const unsigned int  n,
const unsigned int  n_local,
const ParallelType  type = AUTOMATIC 
) [inline]

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

Definition at line 639 of file trilinos_epetra_vector.h.

References EpetraVector< T >::init().

00642 : _destroy_vec_on_exit(true),
00643   myFirstID_(0),
00644   myNumIDs_(0),
00645   myCoefs_(NULL),
00646   nonlocalIDs_(NULL),
00647   nonlocalElementSize_(NULL),
00648   numNonlocalIDs_(0),
00649   allocatedNonlocalLength_(0),
00650   nonlocalCoefs_(NULL),
00651   ignoreNonLocalEntries_(false)
00652 {
00653   this->init(n, n_local, false, type);
00654 }

template<typename T >
EpetraVector< T >::EpetraVector ( const unsigned int  N,
const unsigned int  n_local,
const std::vector< unsigned int > &  ghost,
const ParallelType  type = 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 692 of file trilinos_epetra_vector.h.

References EpetraVector< T >::init().

00696 : _destroy_vec_on_exit(true),
00697   myFirstID_(0),
00698   myNumIDs_(0),
00699   myCoefs_(NULL),
00700   nonlocalIDs_(NULL),
00701   nonlocalElementSize_(NULL),
00702   numNonlocalIDs_(0),
00703   allocatedNonlocalLength_(0),
00704   nonlocalCoefs_(NULL),
00705   ignoreNonLocalEntries_(false)
00706 {
00707   this->init(n, n_local, ghost, false, type);
00708 }

template<typename T >
EpetraVector< T >::EpetraVector ( Epetra_Vector &  v  )  [inline]

Constructor. Creates a EpetraVector assuming you already have a valid Epetra Vec object. In this case, v is NOT destroyed by the EpetraVector 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 EpetraVector.

Definition at line 661 of file trilinos_epetra_vector.h.

References NumericVector< T >::_is_initialized, NumericVector< T >::_type, EpetraVector< T >::_vec, EpetraVector< T >::myCoefs_, EpetraVector< T >::myFirstID_, EpetraVector< T >::myNumIDs_, and libMeshEnums::PARALLEL.

00662 : _destroy_vec_on_exit(false),
00663   myFirstID_(0),
00664   myNumIDs_(0),
00665   myCoefs_(NULL),
00666   nonlocalIDs_(NULL),
00667   nonlocalElementSize_(NULL),
00668   numNonlocalIDs_(0),
00669   allocatedNonlocalLength_(0),
00670   nonlocalCoefs_(NULL),
00671   ignoreNonLocalEntries_(false)
00672 {
00673   _vec = &v;
00674 
00675   this->_type = PARALLEL; // FIXME - need to determine this from v!
00676 
00677   myFirstID_ = _vec->Map().MinMyGID();
00678   myNumIDs_ = _vec->Map().NumMyElements();
00679 
00680   //Currently we impose the restriction that NumVectors==1, so we won't
00681   //need the LDA argument when calling ExtractView. Hence the "dummy" arg.
00682   int dummy;
00683   _vec->ExtractView(&myCoefs_, &dummy);
00684 
00685   this->_is_initialized = true;
00686 }

template<typename T >
EpetraVector< T >::~EpetraVector (  )  [inline]

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

Definition at line 725 of file trilinos_epetra_vector.h.

References EpetraVector< T >::clear().

00726 {
00727   this->clear ();
00728 }


Member Function Documentation

template<typename T >
void EpetraVector< T >::abs (  )  [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 308 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec.

00309 {
00310   _vec->Abs(*_vec);
00311 }

template<typename T >
void EpetraVector< T >::add ( const T  a,
const NumericVector< T > &  v 
) [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 237 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, and EpetraVector< T >::size().

00238 {
00239   const EpetraVector<T>* v = libmesh_cast_ptr<const EpetraVector<T>*>(&v_in);
00240 
00241   libmesh_assert(this->size() == v->size());
00242 
00243   _vec->Update(a_in,*v->_vec, 1.);
00244 }

template<typename T >
void EpetraVector< T >::add ( const NumericVector< T > &  V  )  [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 230 of file trilinos_epetra_vector.C.

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

template<typename T >
void EpetraVector< T >::add ( const T  s  )  [inline, 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 216 of file trilinos_epetra_vector.C.

References NumericVector< T >::_is_closed, and EpetraVector< T >::_vec.

00217 {
00218   const unsigned int nl = _vec->MyLength();
00219 
00220   T * values = _vec->Values();
00221   
00222   for (unsigned int i=0; i<nl; i++)
00223     values[i]+=v_in;
00224 
00225   this->_is_closed = false;
00226 }

template<typename T >
void EpetraVector< T >::add ( const unsigned int  i,
const T  value 
) [inline, virtual]

v(i) += value

Implements NumericVector< T >.

Definition at line 141 of file trilinos_epetra_vector.C.

References NumericVector< T >::_is_closed, EpetraVector< T >::size(), and EpetraVector< T >::SumIntoGlobalValues().

00142 {
00143   int i = static_cast<int> (i_in);
00144   T value = value_in;
00145 
00146   libmesh_assert(i_in<this->size());
00147   
00148   SumIntoGlobalValues(1, &i, &value);
00149 
00150   this->_is_closed = false;
00151 }

template<typename T >
void NumericVector< T >::add_vector ( const NumericVector< T > &  v,
const ShellMatrix< T > &  a 
) [inline, 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().

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

template<typename T >
void EpetraVector< T >::add_vector ( const DenseVector< T > &  V,
const std::vector< unsigned int > &  dof_indices 
) [inline, 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 203 of file trilinos_epetra_vector.C.

References DenseVector< T >::size(), and EpetraVector< T >::SumIntoGlobalValues().

00205 {
00206   libmesh_assert (V_in.size() == dof_indices.size());
00207 
00208   SumIntoGlobalValues(dof_indices.size(),
00209                       (int *)&dof_indices[0],
00210                       &const_cast<DenseVector<T> *>(&V_in)->get_values()[0]);
00211 }

template<typename T >
void EpetraVector< T >::add_vector ( const NumericVector< T > &  V,
const SparseMatrix< T > &  A 
) [inline, 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 182 of file trilinos_epetra_vector.C.

00184 {
00185   libmesh_not_implemented();
00186 
00187 //   const EpetraVector<T>* V = libmesh_cast_ptr<const EpetraVector<T>*>(&V_in);
00188 //   const EpetraMatrix<T>* A = libmesh_cast_ptr<const EpetraMatrix<T>*>(&A_in);
00189 
00190 //   int ierr=0;
00191 
00192 //   A->close();
00193 
00194 //   // The const_cast<> is not elegant, but it is required since Epetra
00195 //   // is not const-correct.  
00196 //   ierr = MatMultAdd(const_cast<EpetraMatrix<T>*>(A)->mat(), V->_vec, _vec, _vec);
00197 //          CHKERRABORT(libMesh::COMM_WORLD,ierr); 
00198 }

template<typename T >
void EpetraVector< T >::add_vector ( const NumericVector< T > &  V,
const std::vector< unsigned int > &  dof_indices 
) [inline, 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 169 of file trilinos_epetra_vector.C.

References NumericVector< T >::size().

00171 {
00172   libmesh_assert (V.size() == dof_indices.size());
00173 
00174   for (unsigned int i=0; i<V.size(); i++)
00175     this->add (dof_indices[i], V(i));
00176 }

template<typename T >
void EpetraVector< T >::add_vector ( const std::vector< T > &  v,
const std::vector< unsigned int > &  dof_indices 
) [inline, 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 156 of file trilinos_epetra_vector.C.

References EpetraVector< T >::SumIntoGlobalValues().

00158 {
00159   libmesh_assert (v.size() == dof_indices.size());
00160 
00161   SumIntoGlobalValues (v.size(),
00162                        (int*) &dof_indices[0],
00163                        const_cast<T*>(&v[0]));
00164 }

template<typename T >
AutoPtr< NumericVector< T > > NumericVector< T >::build ( const SolverPackage  solver_package = libMesh::default_solver_package()  )  [inline, 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().

00042 {
00043   // Build the appropriate vector
00044   switch (solver_package)
00045     {
00046 
00047 
00048 #ifdef LIBMESH_HAVE_LASPACK
00049     case LASPACK_SOLVERS:
00050       {
00051         AutoPtr<NumericVector<T> > ap(new LaspackVector<T>);
00052         return ap;
00053       }
00054 #endif
00055 
00056 
00057 #ifdef LIBMESH_HAVE_PETSC
00058     case PETSC_SOLVERS:
00059       {
00060         AutoPtr<NumericVector<T> > ap(new PetscVector<T>);
00061         return ap;
00062       }
00063 #endif
00064 
00065 
00066 #ifdef LIBMESH_HAVE_TRILINOS
00067     case TRILINOS_SOLVERS:
00068       {
00069         AutoPtr<NumericVector<T> > ap(new EpetraVector<T>);
00070         return ap;
00071       }
00072 #endif
00073 
00074 
00075     default:
00076       AutoPtr<NumericVector<T> > ap(new DistributedVector<T>);
00077       return ap;
00078     }
00079     
00080   AutoPtr<NumericVector<T> > ap(NULL);
00081   return ap;    
00082 }

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

Returns:
the EpetraVector<T> to a pristine state.

Reimplemented from NumericVector< T >.

Definition at line 814 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_destroy_vec_on_exit, NumericVector< T >::_is_closed, NumericVector< T >::_is_initialized, EpetraVector< T >::_vec, and NumericVector< T >::initialized().

Referenced by EpetraVector< T >::~EpetraVector().

00815 {
00816   if ((this->initialized()) && (this->_destroy_vec_on_exit))
00817     delete _vec;
00818 
00819   this->_is_closed = this->_is_initialized = false;
00820 }

template<typename T >
AutoPtr< NumericVector< T > > EpetraVector< T >::clone (  )  const [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 850 of file trilinos_epetra_vector.h.

00851 {
00852   AutoPtr<NumericVector<T> > cloned_vector (new EpetraVector<T>);
00853 
00854   cloned_vector->init(*this, true);
00855 
00856   *cloned_vector = *this;
00857 
00858   return cloned_vector;
00859 }

template<typename T >
void EpetraVector< T >::close (  )  [inline, virtual]

Call the assemble functions

Implements NumericVector< T >.

Definition at line 801 of file trilinos_epetra_vector.h.

References NumericVector< T >::_is_closed, EpetraVector< T >::GlobalAssemble(), and NumericVector< T >::initialized().

Referenced by Problem_Interface::computeF().

00802 {
00803   libmesh_assert (this->initialized());
00804 
00805   this->GlobalAssemble();
00806 
00807   this->_is_closed = true;
00808 }

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

Definition at line 167 of file numeric_vector.C.

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

00169 {
00170   libmesh_assert (this->initialized());
00171   libmesh_assert (other_vector.initialized());
00172   libmesh_assert (this->first_local_index() == other_vector.first_local_index());
00173   libmesh_assert (this->last_local_index()  == other_vector.last_local_index());
00174 
00175   int rvalue     = -1;
00176   unsigned int i = first_local_index();
00177 
00178   do
00179     {
00180       if (( std::abs( (*this)(i).real() - other_vector(i).real() ) > threshold ) ||
00181           ( std::abs( (*this)(i).imag() - other_vector(i).imag() ) > threshold ))
00182         rvalue = i;
00183       else
00184         i++;
00185     }
00186   while (rvalue==-1 && i<this->last_local_index());
00187 
00188   return rvalue;
00189 }

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

Definition at line 140 of file numeric_vector.C.

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

00142 {
00143   libmesh_assert (this->initialized());
00144   libmesh_assert (other_vector.initialized());
00145   libmesh_assert (this->first_local_index() == other_vector.first_local_index());
00146   libmesh_assert (this->last_local_index()  == other_vector.last_local_index());
00147 
00148   int rvalue     = -1;
00149   unsigned int i = first_local_index();
00150 
00151   do
00152     {
00153       if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
00154         rvalue = i;
00155       else
00156         i++;
00157     }
00158   while (rvalue==-1 && i<last_local_index());
00159 
00160   return rvalue;
00161 }

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

Definition at line 114 of file numeric_vector.C.

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

00116 {
00117   libmesh_assert (this->initialized());
00118   libmesh_assert (other_vector.initialized());
00119   libmesh_assert (this->first_local_index() == other_vector.first_local_index());
00120   libmesh_assert (this->last_local_index()  == other_vector.last_local_index());
00121 
00122   int rvalue     = -1;
00123   unsigned int i = first_local_index();
00124 
00125   do
00126     {
00127       if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
00128         rvalue = i;
00129       else
00130         i++;
00131     }
00132   while (rvalue==-1 && i<last_local_index());
00133 
00134   return rvalue;
00135 }

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

Definition at line 89 of file numeric_vector.C.

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

00091 {
00092   libmesh_assert (this->initialized());
00093   libmesh_assert (other_vector.initialized());
00094   libmesh_assert (this->first_local_index() == other_vector.first_local_index());
00095   libmesh_assert (this->last_local_index()  == other_vector.last_local_index());
00096 
00097   int rvalue     = -1;
00098   unsigned int i = first_local_index();
00099 
00100   do
00101     {
00102       if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
00103         rvalue = i;
00104       else
00105         i++;
00106     }
00107   while (rvalue==-1 && i<last_local_index());
00108 
00109   return rvalue;
00110 }

template<typename T >
virtual int NumericVector< T >::compare ( const NumericVector< T > &  other_vector,
const Real  threshold = 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<typename T >
void EpetraVector< T >::create_subvector ( NumericVector< T > &  subvector,
const std::vector< unsigned int > &  rows 
) const [inline, virtual]

Creates a "subvector" from this vector using the rows indices of the "rows" array.

Reimplemented from NumericVector< T >.

Definition at line 644 of file trilinos_epetra_vector.C.

00646 {
00647   libmesh_not_implemented();
00648 
00649 //   // Epetra data structures
00650 //   IS parent_is, subvector_is;
00651 //   VecScatter scatter;
00652 //   int ierr = 0;
00653   
00654 //   // Make sure the passed int subvector is really a EpetraVector
00655 //   EpetraVector<T>* epetra_subvector = libmesh_cast_ptr<EpetraVector<T>*>(&subvector);
00656 //   libmesh_assert(epetra_subvector != NULL);
00657   
00658 //   // If the epetra_subvector is already initialized, we assume that the
00659 //   // user has already allocated the *correct* amount of space for it.
00660 //   // If not, we use the appropriate Epetra routines to initialize it.
00661 //   if (!epetra_subvector->initialized())
00662 //     {
00663 //       // Initialize the epetra_subvector to have enough space to hold
00664 //       // the entries which will be scattered into it.  Note: such an
00665 //       // init() function (where we let Epetra decide the number of local
00666 //       // entries) is not currently offered by the EpetraVector
00667 //       // class.  Should we differentiate here between sequential and
00668 //       // parallel vector creation based on libMesh::n_processors() ?
00669 //       ierr = VecCreateMPI(libMesh::COMM_WORLD,
00670 //                        EPETRA_DECIDE,          // n_local
00671 //                        rows.size(),           // n_global
00672 //                        &(epetra_subvector->_vec)); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00673 
00674 //       ierr = VecSetFromOptions (epetra_subvector->_vec); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00675 
00676 //       // Mark the subvector as initialized
00677 //       epetra_subvector->_is_initialized = true;
00678 //     }
00679   
00680 //   // Use iota to fill an array with entries [0,1,2,3,4,...rows.size()]
00681 //   std::vector<int> idx(rows.size());
00682 //   Utility::iota (idx.begin(), idx.end(), 0);
00683 
00684 //   // Construct index sets
00685 //   ierr = ISCreateGeneral(libMesh::COMM_WORLD,
00686 //                       rows.size(),
00687 //                       (int*) &rows[0],
00688 //                       &parent_is); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00689 
00690 //   ierr = ISCreateGeneral(libMesh::COMM_WORLD,
00691 //                       rows.size(),
00692 //                       (int*) &idx[0],
00693 //                       &subvector_is); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00694 
00695 //   // Construct the scatter object
00696 //   ierr = VecScatterCreate(this->_vec,
00697 //                        parent_is,
00698 //                        epetra_subvector->_vec,
00699 //                        subvector_is,
00700 //                        &scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00701 
00702 //   // Actually perform the scatter
00703 // #if EPETRA_VERSION_LESS_THAN(2,3,3)
00704 //   ierr = VecScatterBegin(this->_vec,
00705 //                       epetra_subvector->_vec,
00706 //                       INSERT_VALUES,
00707 //                       SCATTER_FORWARD,
00708 //                       scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00709 
00710 //   ierr = VecScatterEnd(this->_vec,
00711 //                     epetra_subvector->_vec,
00712 //                     INSERT_VALUES,
00713 //                     SCATTER_FORWARD,
00714 //                     scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00715 // #else
00716 //   // API argument order change in Epetra 2.3.3
00717 //   ierr = VecScatterBegin(scatter,
00718 //                       this->_vec,
00719 //                       epetra_subvector->_vec,
00720 //                       INSERT_VALUES,
00721 //                       SCATTER_FORWARD); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00722 
00723 //   ierr = VecScatterEnd(scatter,
00724 //                     this->_vec,
00725 //                     epetra_subvector->_vec,
00726 //                     INSERT_VALUES,
00727 //                     SCATTER_FORWARD); CHKERRABORT(libMesh::COMM_WORLD,ierr);
00728 // #endif
00729   
00730 //   // Clean up 
00731 //   ierr = ISDestroy(parent_is);       CHKERRABORT(libMesh::COMM_WORLD,ierr);
00732 //   ierr = ISDestroy(subvector_is);    CHKERRABORT(libMesh::COMM_WORLD,ierr);
00733 //   ierr = VecScatterDestroy(scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr); 
00734 }

template<typename T >
void EpetraVector< T >::destroyNonlocalData (  )  [inline, private]

Definition at line 1040 of file trilinos_epetra_vector.C.

References EpetraVector< T >::allocatedNonlocalLength_, EpetraVector< T >::nonlocalCoefs_, EpetraVector< T >::nonlocalElementSize_, EpetraVector< T >::nonlocalIDs_, and EpetraVector< T >::numNonlocalIDs_.

Referenced by EpetraVector< T >::FEoperatorequals(), and EpetraVector< T >::GlobalAssemble().

01041 {
01042   if (allocatedNonlocalLength_ > 0) {
01043     delete [] nonlocalIDs_;
01044     delete [] nonlocalElementSize_;
01045     nonlocalIDs_ = NULL;
01046     nonlocalElementSize_ = NULL;
01047     for(int i=0; i<numNonlocalIDs_; ++i) {
01048       delete [] nonlocalCoefs_[i];
01049     }
01050     delete [] nonlocalCoefs_;
01051     nonlocalCoefs_ = NULL;
01052     numNonlocalIDs_ = 0;
01053     allocatedNonlocalLength_ = 0;
01054   }
01055   return;
01056 }

template<typename T >
T EpetraVector< T >::dot ( const NumericVector< T > &  V  )  const [inline, virtual]

Computes the dot product, p = U.V

Implements NumericVector< T >.

Definition at line 315 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec.

00316 {
00317   const EpetraVector<T>* V = libmesh_cast_ptr<const EpetraVector<T>*>(&V_in);
00318 
00319   T result=0.0;
00320 
00321   _vec->Dot(*V->_vec, &result);
00322 
00323   return result;
00324 }

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

Returns:
the element U(i)

Definition at line 322 of file numeric_vector.h.

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

template<typename T >
void EpetraVector< T >::FEoperatorequals ( const EpetraVector< T > &  source  )  [inline, private]

Definition at line 1013 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, EpetraVector< T >::allocatedNonlocalLength_, EpetraVector< T >::destroyNonlocalData(), EpetraVector< T >::nonlocalCoefs_, EpetraVector< T >::nonlocalElementSize_, EpetraVector< T >::nonlocalIDs_, and EpetraVector< T >::numNonlocalIDs_.

01014 {
01015   (*_vec) = *(source._vec);
01016 
01017   destroyNonlocalData();
01018 
01019   if (source.allocatedNonlocalLength_ > 0) {
01020     allocatedNonlocalLength_ = source.allocatedNonlocalLength_;
01021     numNonlocalIDs_ = source.numNonlocalIDs_;
01022     nonlocalIDs_ = new int[allocatedNonlocalLength_];
01023     nonlocalElementSize_ = new int[allocatedNonlocalLength_];
01024     nonlocalCoefs_ = new double*[allocatedNonlocalLength_];
01025     for(int i=0; i<numNonlocalIDs_; ++i) {
01026       int elemSize = source.nonlocalElementSize_[i];
01027       nonlocalCoefs_[i] = new double[elemSize];
01028       nonlocalIDs_[i] = source.nonlocalIDs_[i];
01029       nonlocalElementSize_[i] = elemSize;
01030       for(int j=0; j<elemSize; ++j) {
01031         nonlocalCoefs_[i][j] = source.nonlocalCoefs_[i][j];
01032       }
01033     }
01034   }
01035 }

template<typename T >
unsigned int EpetraVector< 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 885 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec, and NumericVector< T >::initialized().

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

00886 {
00887   libmesh_assert (this->initialized());
00888   
00889   return _vec->Map().MinMyGID();
00890 }

template<typename T >
void NumericVector< T >::get ( const std::vector< unsigned int > &  index,
std::vector< T > &  values 
) const [inline, virtual, inherited]

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

Reimplemented in PetscVector< T >.

Definition at line 722 of file numeric_vector.h.

00723 {
00724   const unsigned int num = index.size();
00725   values.resize(num);
00726   for(unsigned int i=0; i<num; i++)
00727     {
00728       values[i] = (*this)(index[i]);
00729     }
00730 }

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 QuadratureRules::name().

Referenced by ReferenceCounter::print_info().

00046 {
00047 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
00048 
00049   std::ostringstream out;
00050   
00051   out << '\n'
00052       << " ---------------------------------------------------------------------------- \n"
00053       << "| Reference count information                                                |\n"
00054       << " ---------------------------------------------------------------------------- \n";
00055   
00056   for (Counts::iterator it = _counts.begin();
00057        it != _counts.end(); ++it)
00058     {
00059       const std::string name(it->first);
00060       const unsigned int creations    = it->second.first;
00061       const unsigned int destructions = it->second.second;
00062 
00063       out << "| " << name << " reference count information:\n"
00064           << "|  Creations:    " << creations    << '\n'
00065           << "|  Destructions: " << destructions << '\n';
00066     }
00067   
00068   out << " ---------------------------------------------------------------------------- \n";
00069 
00070   return out.str();
00071 
00072 #else
00073 
00074   return "";
00075   
00076 #endif
00077 }

template<typename T >
int EpetraVector< T >::GlobalAssemble ( Epetra_CombineMode  mode = Add  )  [inline, private]

Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this vector at construction time. Data imported from other processors is stored on the owning processor with a "sumInto" or accumulate operation. This is a collective method -- every processor must enter it before any will complete it.

Definition at line 968 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, EpetraVector< T >::destroyNonlocalData(), EpetraVector< T >::ignoreNonLocalEntries_, EpetraVector< T >::nonlocalCoefs_, EpetraVector< T >::nonlocalElementSize_, EpetraVector< T >::nonlocalIDs_, and EpetraVector< T >::numNonlocalIDs_.

Referenced by EpetraVector< T >::close().

00969 {
00970   //In this method we need to gather all the non-local (overlapping) data
00971   //that's been input on each processor, into the (probably) non-overlapping
00972   //distribution defined by the map that 'this' vector was constructed with.
00973 
00974   //We don't need to do anything if there's only one processor or if
00975   //ignoreNonLocalEntries_ is true.
00976   if (_vec->Map().Comm().NumProc() < 2 || ignoreNonLocalEntries_) {
00977     return(0);
00978   }
00979 
00980 
00981   
00982   //First build a map that describes the data in nonlocalIDs_/nonlocalCoefs_.
00983   //We'll use the arbitrary distribution constructor of Map.
00984 
00985   Epetra_BlockMap sourceMap(-1, numNonlocalIDs_,
00986                             nonlocalIDs_, nonlocalElementSize_,
00987                             _vec->Map().IndexBase(), _vec->Map().Comm());
00988 
00989   //Now build a vector to hold our nonlocalCoefs_, and to act as the source-
00990   //vector for our import operation.
00991   Epetra_MultiVector nonlocalVector(sourceMap, 1);
00992 
00993   int i,j;
00994   for(i=0; i<numNonlocalIDs_; ++i) {
00995     for(j=0; j<nonlocalElementSize_[i]; ++j) {
00996       nonlocalVector.ReplaceGlobalValue(nonlocalIDs_[i], j, 0,
00997                                         nonlocalCoefs_[i][j]);
00998     }
00999   }
01000 
01001   Epetra_Export exporter(sourceMap, _vec->Map());
01002 
01003   EPETRA_CHK_ERR( _vec->Export(nonlocalVector, exporter, mode) );
01004 
01005   destroyNonlocalData();
01006 
01007   return(0);
01008 }

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, and Threads::spin_mtx.

Referenced by ReferenceCountedObject< SparseMatrix< T > >::ReferenceCountedObject().

00150 {
00151   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00152   std::pair<unsigned int, unsigned int>& p = _counts[name];
00153 
00154   p.first++;
00155 }

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, and Threads::spin_mtx.

Referenced by ReferenceCountedObject< SparseMatrix< T > >::~ReferenceCountedObject().

00168 {
00169   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00170   std::pair<unsigned int, unsigned int>& p = _counts[name];
00171 
00172   p.second++;
00173 }

template<class T >
void EpetraVector< T >::init ( const NumericVector< T > &  other,
const bool  fast = 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 715 of file trilinos_epetra_vector.h.

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

00717 {
00718   this->init(other.size(),other.local_size(),fast,other.type());
00719 }

template<typename T >
void EpetraVector< T >::init ( const unsigned int  n,
const unsigned int  n_local,
const std::vector< unsigned int > &  ,
const bool  fast = false,
const ParallelType  type = AUTOMATIC 
) [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 776 of file trilinos_epetra_vector.h.

References EpetraVector< T >::init().

00781 {
00782   // TODO: we shouldn't ignore the ghost sparsity pattern
00783   this->init(n, n_local, fast, type);
00784 }

template<typename T >
void EpetraVector< T >::init ( const unsigned int  N,
const bool  fast = false,
const ParallelType  type = AUTOMATIC 
) [inline, virtual]

call init with n_local = N,

Implements NumericVector< T >.

Definition at line 790 of file trilinos_epetra_vector.h.

References EpetraVector< T >::init().

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

template<typename T >
void EpetraVector< T >::init ( const unsigned int  N,
const unsigned int  n_local,
const bool  fast = false,
const ParallelType  type = 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 734 of file trilinos_epetra_vector.h.

References NumericVector< T >::_is_initialized, EpetraVector< T >::_map, NumericVector< T >::_type, EpetraVector< T >::_vec, libMeshEnums::AUTOMATIC, libMesh::COMM_WORLD, EpetraVector< T >::myCoefs_, EpetraVector< T >::myFirstID_, EpetraVector< T >::myNumIDs_, libMeshEnums::PARALLEL, libMeshEnums::SERIAL, and EpetraVector< T >::zero().

Referenced by EpetraVector< T >::EpetraVector(), and EpetraVector< T >::init().

00738 {
00739   if (type == AUTOMATIC)
00740     {
00741       if (n == n_local)
00742         this->_type = SERIAL;
00743       else
00744         this->_type = PARALLEL;
00745     }
00746   else
00747     this->_type = type;
00748 
00749   libmesh_assert ((this->_type==SERIAL && n==n_local) ||
00750                   this->_type==PARALLEL);
00751 
00752   _map = new Epetra_Map(n, 
00753                         n_local, 
00754                         0, 
00755                         Epetra_MpiComm (libMesh::COMM_WORLD));
00756               
00757   _vec = new Epetra_Vector(*_map);
00758 
00759   myFirstID_ = _vec->Map().MinMyGID();
00760   myNumIDs_ = _vec->Map().NumMyElements();
00761 
00762   //Currently we impose the restriction that NumVectors==1, so we won't
00763   //need the LDA argument when calling ExtractView. Hence the "dummy" arg.
00764   int dummy;
00765   _vec->ExtractView(&myCoefs_, &dummy);
00766   
00767   this->_is_initialized = true;
00768   
00769   if (fast == false)
00770     this->zero ();
00771 }

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 PetscVector< T >::_get_array(), PetscVector< T >::_restore_array(), LaspackVector< T >::abs(), DistributedVector< T >::abs(), LaspackVector< T >::add(), DistributedVector< T >::add(), DistributedVector< T >::add_vector(), EpetraVector< T >::clear(), PetscVector< T >::clear(), LaspackVector< T >::clear(), EpetraVector< T >::close(), LaspackVector< T >::close(), DistributedVector< T >::close(), NumericVector< T >::compare(), PetscVector< T >::create_subvector(), LaspackVector< T >::dot(), EpetraVector< T >::first_local_index(), PetscVector< T >::first_local_index(), LaspackVector< T >::first_local_index(), DistributedVector< T >::first_local_index(), PetscVector< T >::init(), LaspackVector< T >::init(), DistributedVector< T >::init(), DistributedVector< T >::insert(), DistributedVector< T >::l1_norm(), DistributedVector< T >::l2_norm(), EpetraVector< T >::last_local_index(), PetscVector< T >::last_local_index(), LaspackVector< T >::last_local_index(), DistributedVector< T >::last_local_index(), DistributedVector< T >::linfty_norm(), EpetraVector< T >::local_size(), PetscVector< T >::local_size(), LaspackVector< T >::local_size(), DistributedVector< T >::local_size(), DistributedVector< T >::localize(), DistributedVector< T >::localize_to_one(), PetscVector< T >::map_global_to_local_index(), EpetraVector< T >::max(), LaspackVector< T >::max(), DistributedVector< T >::max(), EpetraVector< T >::min(), LaspackVector< T >::min(), DistributedVector< T >::min(), EpetraVector< T >::operator()(), LaspackVector< T >::operator()(), DistributedVector< T >::operator()(), DistributedVector< T >::operator+=(), DistributedVector< T >::operator-=(), LaspackVector< T >::operator=(), DistributedVector< T >::operator=(), NumericVector< T >::print(), NumericVector< T >::print_global(), LaspackVector< T >::scale(), DistributedVector< T >::scale(), LaspackVector< T >::set(), DistributedVector< T >::set(), EpetraVector< T >::size(), PetscVector< T >::size(), LaspackVector< T >::size(), DistributedVector< T >::size(), DistributedVector< T >::sum(), EpetraVector< T >::zero(), LaspackVector< T >::zero(), and DistributedVector< T >::zero().

00109 { return _is_initialized; }

template<typename T >
int EpetraVector< T >::inputNonlocalValue ( int  GID,
double  value,
bool  accumulate 
) [inline, private]

Definition at line 868 of file trilinos_epetra_vector.C.

References EpetraVector< T >::allocatedNonlocalLength_, EpetraVector< T >::nonlocalCoefs_, EpetraVector< T >::nonlocalElementSize_, EpetraVector< T >::nonlocalIDs_, and EpetraVector< T >::numNonlocalIDs_.

Referenced by EpetraVector< T >::inputValues().

00869 {
00870   int insertPoint = -1;
00871 
00872   //find offset of GID in nonlocalIDs_
00873   int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_,
00874                                          insertPoint);
00875   if (offset >= 0) {
00876     //if offset >= 0
00877     //  put value in nonlocalCoefs_[offset][0]
00878 
00879     if (accumulate) {
00880       nonlocalCoefs_[offset][0] += value;
00881     }
00882     else {
00883       nonlocalCoefs_[offset][0] = value;
00884     }
00885   }
00886   else {
00887     //else
00888     //  insert GID in nonlocalIDs_
00889     //  insert 1   in nonlocalElementSize_
00890     //  insert value in nonlocalCoefs_
00891 
00892     int tmp1 = numNonlocalIDs_;
00893     int tmp2 = allocatedNonlocalLength_;
00894     int tmp3 = allocatedNonlocalLength_;
00895     EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_,
00896                                        tmp1, tmp2) );
00897     --tmp1;
00898     EPETRA_CHK_ERR( Epetra_Util_insert(1, insertPoint, nonlocalElementSize_,
00899                                        tmp1, tmp3) );
00900     double* values = new double[1];
00901     values[0] = value;
00902     EPETRA_CHK_ERR( Epetra_Util_insert(values, insertPoint, nonlocalCoefs_,
00903                                        numNonlocalIDs_, allocatedNonlocalLength_) );
00904   }
00905 
00906   return(0);
00907 }

template<typename T >
int EpetraVector< T >::inputNonlocalValues ( int  GID,
int  numValues,
const double *  values,
bool  accumulate 
) [inline, private]

Definition at line 911 of file trilinos_epetra_vector.C.

References EpetraVector< T >::allocatedNonlocalLength_, EpetraVector< T >::nonlocalCoefs_, EpetraVector< T >::nonlocalElementSize_, EpetraVector< T >::nonlocalIDs_, and EpetraVector< T >::numNonlocalIDs_.

Referenced by EpetraVector< T >::inputValues().

00913 {
00914   int insertPoint = -1;
00915 
00916   //find offset of GID in nonlocalIDs_
00917   int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_,
00918                                          insertPoint);
00919   if (offset >= 0) {
00920     //if offset >= 0
00921     //  put value in nonlocalCoefs_[offset][0]
00922 
00923     if (numValues != nonlocalElementSize_[offset]) {
00924       cerr << "Epetra_FEVector ERROR: block-size for GID " << GID << " is "
00925            << numValues<<" which doesn't match previously set block-size of "
00926            << nonlocalElementSize_[offset] << endl;
00927       return(-1);
00928     }
00929 
00930     if (accumulate) {
00931       for(int j=0; j<numValues; ++j) {
00932         nonlocalCoefs_[offset][j] += values[j];
00933       }
00934     }
00935     else {
00936       for(int j=0; j<numValues; ++j) {
00937         nonlocalCoefs_[offset][j] = values[j];
00938       }
00939     }
00940   }
00941   else {
00942     //else
00943     //  insert GID in nonlocalIDs_
00944     //  insert numValues   in nonlocalElementSize_
00945     //  insert values in nonlocalCoefs_
00946 
00947     int tmp1 = numNonlocalIDs_;
00948     int tmp2 = allocatedNonlocalLength_;
00949     int tmp3 = allocatedNonlocalLength_;
00950     EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_,
00951                                        tmp1, tmp2) );
00952     --tmp1;
00953     EPETRA_CHK_ERR( Epetra_Util_insert(numValues, insertPoint, nonlocalElementSize_,
00954                                        tmp1, tmp3) );
00955     double* newvalues = new double[numValues];
00956     for(int j=0; j<numValues; ++j) {
00957       newvalues[j] = values[j];
00958     }
00959     EPETRA_CHK_ERR( Epetra_Util_insert(newvalues, insertPoint, nonlocalCoefs_,
00960                                        numNonlocalIDs_, allocatedNonlocalLength_) );
00961   }
00962 
00963   return(0);
00964 }

template<typename T >
int EpetraVector< T >::inputValues ( int  numIDs,
const int *  GIDs,
const int *  numValuesPerID,
const double *  values,
bool  accumulate 
) [inline, private]

Definition at line 833 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, EpetraVector< T >::ignoreNonLocalEntries_, and EpetraVector< T >::inputNonlocalValues().

00838 {
00839   int offset=0;
00840   for(int i=0; i<numIDs; ++i) {
00841     int numValues = numValuesPerID[i];
00842     if (_vec->Map().MyGID(GIDs[i])) {
00843       if (accumulate) {
00844         for(int j=0; j<numValues; ++j) {
00845           _vec->SumIntoGlobalValue(GIDs[i], j, 0, values[offset+j]);
00846         }
00847       }
00848       else {
00849         for(int j=0; j<numValues; ++j) {
00850           _vec->ReplaceGlobalValue(GIDs[i], j, 0, values[offset+j]);
00851         }
00852       }
00853     }
00854     else {
00855       if (!ignoreNonLocalEntries_) {
00856         EPETRA_CHK_ERR( inputNonlocalValues(GIDs[i], numValues,
00857                                             &(values[offset]), accumulate) );
00858       }
00859     }
00860     offset += numValues;
00861   }
00862 
00863   return(0);
00864 }

template<typename T >
int EpetraVector< T >::inputValues ( int  numIDs,
const int *  GIDs,
const double *  values,
bool  accumulate 
) [inline, private]

Definition at line 804 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, EpetraVector< T >::ignoreNonLocalEntries_, and EpetraVector< T >::inputNonlocalValue().

Referenced by EpetraVector< T >::ReplaceGlobalValues(), and EpetraVector< T >::SumIntoGlobalValues().

00808 {
00809  //Important note!! This method assumes that there is only 1 point
00810  //associated with each element.
00811 
00812   for(int i=0; i<numIDs; ++i) {
00813     if (_vec->Map().MyGID(GIDs[i])) {
00814       if (accumulate) {
00815         _vec->SumIntoGlobalValue(GIDs[i], 0, 0, values[i]);
00816       }
00817       else {
00818         _vec->ReplaceGlobalValue(GIDs[i], 0, 0, values[i]);
00819       }
00820     }
00821     else {
00822       if (!ignoreNonLocalEntries_) {
00823         EPETRA_CHK_ERR( inputNonlocalValue(GIDs[i], values[i], accumulate) );
00824       }
00825     }
00826   }
00827 
00828   return(0);
00829 }

template<typename T >
void EpetraVector< T >::insert ( const DenseSubVector< T > &  V,
const std::vector< unsigned int > &  dof_indices 
) [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 290 of file trilinos_epetra_vector.C.

References DenseVectorBase< T >::size().

00292 {
00293   libmesh_assert (v.size() == dof_indices.size());
00294   
00295   for (unsigned int i=0; i < v.size(); ++i)
00296     this->set (dof_indices[i], v(i));
00297 }

template<typename T >
void EpetraVector< T >::insert ( const DenseVector< T > &  V,
const std::vector< unsigned int > &  dof_indices 
) [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 275 of file trilinos_epetra_vector.C.

References EpetraVector< T >::ReplaceGlobalValues(), and DenseVector< T >::size().

00277 {
00278   libmesh_assert (v.size() == dof_indices.size());
00279   
00280   std::vector<T> &vals = const_cast<DenseVector<T>&>(v).get_values();
00281   
00282   ReplaceGlobalValues (v.size(),
00283                        (int*) &dof_indices[0],
00284                        &vals[0]);
00285 }

template<typename T >
void EpetraVector< T >::insert ( const NumericVector< T > &  V,
const std::vector< unsigned int > &  dof_indices 
) [inline, 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 262 of file trilinos_epetra_vector.C.

References NumericVector< T >::size().

00264 {
00265   libmesh_assert (V.size() == dof_indices.size());
00266 
00267   // TODO: If V is an EpetraVector this can be optimized
00268   for (unsigned int i=0; i<V.size(); i++)
00269     this->set (dof_indices[i], V(i));
00270 }

template<typename T >
void EpetraVector< T >::insert ( const std::vector< T > &  v,
const std::vector< unsigned int > &  dof_indices 
) [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 249 of file trilinos_epetra_vector.C.

References EpetraVector< T >::ReplaceGlobalValues().

00251 {
00252   libmesh_assert (v.size() == dof_indices.size());
00253 
00254   ReplaceGlobalValues (v.size(),
00255                        (int*) &dof_indices[0],
00256                        const_cast<T*>(&v[0]));
00257 }

template<typename T >
Real EpetraVector< T >::l1_norm (  )  const [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 64 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, and NumericVector< T >::closed().

00065 {
00066   libmesh_assert(this->closed());
00067 
00068   Real value;
00069 
00070   _vec->Norm1(&value);
00071   
00072   return value;
00073 }

template<typename T >
Real EpetraVector< T >::l2_norm (  )  const [inline, 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 76 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, and NumericVector< T >::closed().

00077 {
00078   libmesh_assert(this->closed());
00079   
00080   Real value;
00081 
00082   _vec->Norm2(&value);
00083   
00084   return value;
00085 }

template<typename T >
unsigned int EpetraVector< 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 896 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec, and NumericVector< T >::initialized().

Referenced by EpetraVector< T >::operator()().

00897 {
00898   libmesh_assert (this->initialized());
00899   
00900   return _vec->Map().MaxMyGID()+1;
00901 }

template<typename T >
Real EpetraVector< T >::linfty_norm (  )  const [inline, 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 88 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, and NumericVector< T >::closed().

00089 {
00090   libmesh_assert(this->closed());
00091   
00092   Real value;
00093 
00094   _vec->NormInf(&value);
00095   
00096   return value;
00097 }

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

Returns:
the local size of the vector (index_stop-index_start)

Implements NumericVector< T >.

Definition at line 876 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec, and NumericVector< T >::initialized().

Referenced by EpetraVector< T >::localize(), EpetraVector< T >::localize_to_one(), and EpetraVector< T >::operator=().

00877 {
00878   libmesh_assert (this->initialized());
00879   
00880   return _vec->MyLength();
00881 }

template<typename T >
void EpetraVector< T >::localize ( const unsigned int  first_local_idx,
const unsigned int  last_local_idx,
const std::vector< unsigned int > &  send_list 
) [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 447 of file trilinos_epetra_vector.C.

00450 {
00451   libmesh_not_implemented();
00452 
00453 //   // Only good for serial vectors.
00454 //   libmesh_assert (this->size() == this->local_size());
00455 //   libmesh_assert (last_local_idx > first_local_idx);
00456 //   libmesh_assert (send_list.size() <= this->size());
00457 //   libmesh_assert (last_local_idx < this->size());
00458   
00459 //   const unsigned int size       = this->size();
00460 //   const unsigned int local_size = (last_local_idx - first_local_idx + 1);
00461 //   int ierr=0;  
00462   
00463 //   // Don't bother for serial cases
00464 //   if ((first_local_idx == 0) &&
00465 //       (local_size == size))
00466 //     return;
00467   
00468   
00469 //   // Build a parallel vector, initialize it with the local
00470 //   // parts of (*this)
00471 //   EpetraVector<T> parallel_vec;
00472 
00473 //   parallel_vec.init (size, local_size, true, PARALLEL);
00474 
00475 
00476 //   // Copy part of *this into the parallel_vec
00477 //   {
00478 //     IS is;
00479 //     VecScatter scatter;
00480 
00481 //     // Create idx, idx[i] = i+first_local_idx;
00482 //     std::vector<int> idx(local_size);
00483 //     Utility::iota (idx.begin(), idx.end(), first_local_idx);
00484 
00485 //     // Create the index set & scatter object
00486 //     ierr = ISCreateGeneral(libMesh::COMM_WORLD, local_size, &idx[0], &is); 
00487 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00488 
00489 //     ierr = VecScatterCreate(_vec,              is,
00490 //                          parallel_vec._vec, is,
00491 //                          &scatter);
00492 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00493 
00494 //     // Perform the scatter
00495 // #if EPETRA_VERSION_LESS_THAN(2,3,3)
00496 
00497 //     ierr = VecScatterBegin(_vec, parallel_vec._vec, INSERT_VALUES,
00498 //                         SCATTER_FORWARD, scatter);
00499 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00500   
00501 //     ierr = VecScatterEnd  (_vec, parallel_vec._vec, INSERT_VALUES,
00502 //                         SCATTER_FORWARD, scatter);
00503 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00504 
00505 // #else
00506            
00507 //       // API argument order change in Epetra 2.3.3
00508 //     ierr = VecScatterBegin(scatter, _vec, parallel_vec._vec,
00509 //                         INSERT_VALUES, SCATTER_FORWARD);
00510 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00511   
00512 //     ierr = VecScatterEnd  (scatter, _vec, parallel_vec._vec,
00513 //                         INSERT_VALUES, SCATTER_FORWARD);
00514 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00515            
00516 // #endif
00517 
00518 //     // Clean up
00519 //     ierr = ISDestroy (is);
00520 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00521   
00522 //     ierr = VecScatterDestroy(scatter);
00523 //            CHKERRABORT(libMesh::COMM_WORLD,ierr);
00524 //   }
00525 
00526 //   // localize like normal
00527 //   parallel_vec.close();
00528 //   parallel_vec.localize (*this, send_list);
00529 //   this->close();
00530 }

template<typename T >
void EpetraVector< T >::localize ( NumericVector< T > &  v_local,
const std::vector< unsigned int > &  send_list 
) const [inline, 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 426 of file trilinos_epetra_vector.C.

References EpetraVector< T >::localize().

00428 {
00429   // TODO: optimize to sync only the send list values
00430   this->localize(v_local_in);
00431   
00432 //   EpetraVector<T>* v_local =
00433 //   libmesh_cast_ptr<EpetraVector<T>*>(&v_local_in);
00434 
00435 //   libmesh_assert (this->_map.get() != NULL);
00436 //   libmesh_assert (v_local->_map.get() != NULL);
00437 //   libmesh_assert (v_local->local_size() == this->size());
00438 //   libmesh_assert (send_list.size() <= v_local->size());
00439 
00440 //   Epetra_Import importer (*v_local->_map, *this->_map);
00441   
00442 //   v_local->_vec->Import (*this->_vec, importer, Insert);
00443 }

template<typename T >
void EpetraVector< T >::localize ( NumericVector< T > &  v_local  )  const [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 412 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_map, and EpetraVector< T >::_vec.

00413 {
00414   EpetraVector<T>* v_local = libmesh_cast_ptr<EpetraVector<T>*>(&v_local_in);
00415 
00416   Epetra_Map rootMap = Epetra_Util::Create_Root_Map( *_map, -1);
00417   v_local->_vec->ReplaceMap(rootMap);
00418 
00419   Epetra_Import importer(v_local->_vec->Map(), *_map);
00420   v_local->_vec->Import(*_vec, importer, Insert);
00421 }

template<typename T >
void EpetraVector< T >::localize ( std::vector< T > &  v_local  )  const [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 535 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, EpetraVector< T >::local_size(), and EpetraVector< T >::size().

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

00536 {
00537   // This function must be run on all processors at once
00538   parallel_only();
00539 
00540   const unsigned int n  = this->size();
00541   const unsigned int nl = this->local_size();
00542 
00543   libmesh_assert (this->_vec != NULL);
00544 
00545   v_local.clear();
00546   v_local.reserve(n);
00547 
00548   // build up my local part
00549   for (unsigned int i=0; i<nl; i++)
00550     v_local.push_back((*this->_vec)[i]);
00551 
00552   Parallel::allgather (v_local);
00553 }

template<typename T >
void EpetraVector< T >::localize_to_one ( std::vector< T > &  v_local,
const unsigned int  proc_id = 0 
) const [inline, 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 >.

Definition at line 558 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, EpetraVector< T >::local_size(), libMesh::n_processors(), and EpetraVector< T >::size().

00560 {
00561   // This function must be run on all processors at once
00562   parallel_only();
00563 
00564   const unsigned int n  = this->size();
00565   const unsigned int nl = this->local_size();
00566 
00567   libmesh_assert (pid < libMesh::n_processors());
00568   libmesh_assert (this->_vec != NULL);
00569   
00570   v_local.clear();
00571   v_local.reserve(n);
00572   
00573   
00574   // build up my local part
00575   for (unsigned int i=0; i<nl; i++)
00576     v_local.push_back((*this->_vec)[i]);
00577 
00578   Parallel::gather (pid, v_local);
00579 }

template<typename T >
Real EpetraVector< 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 934 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec, and NumericVector< T >::initialized().

00935 {
00936   libmesh_assert (this->initialized());
00937 
00938   T value;
00939 
00940   _vec->MaxValue(&value);
00941 
00942   return value;
00943 }

template<typename T >
Real EpetraVector< 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 919 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec, and NumericVector< T >::initialized().

00920 {
00921   libmesh_assert (this->initialized());
00922 
00923   T value;
00924 
00925   _vec->MinValue(&value);
00926 
00927   return value;
00928 }

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.

00077   { return _n_objects; }

template<typename T >
T EpetraVector< T >::operator() ( const unsigned int  i  )  const [inline, virtual]

Access components, returns U(i).

Implements NumericVector< T >.

Definition at line 906 of file trilinos_epetra_vector.h.

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

00907 {
00908   libmesh_assert (this->initialized());
00909   libmesh_assert ( ((i >= this->first_local_index()) &&
00910                     (i <  this->last_local_index())) );
00911 
00912   return (*_vec)[i-this->first_local_index()];
00913 }

template<typename T >
NumericVector<T>& NumericVector< T >::operator*= ( const T  a  )  [inline, inherited]

Multiplication operator. Equivalent to U.scale(a)

Definition at line 348 of file numeric_vector.h.

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

template<typename T >
NumericVector< T > & EpetraVector< T >::operator+= ( const NumericVector< T > &  V  )  [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 101 of file trilinos_epetra_vector.C.

References NumericVector< T >::closed().

00102 {
00103   libmesh_assert(this->closed());
00104   
00105   this->add(1., v);
00106   
00107   return *this;
00108 }

template<typename T >
NumericVector< T > & EpetraVector< T >::operator-= ( const NumericVector< T > &  V  )  [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 114 of file trilinos_epetra_vector.C.

References NumericVector< T >::closed().

00115 {
00116   libmesh_assert(this->closed());
00117   
00118   this->add(-1., v);
00119   
00120   return *this;
00121 }

template<typename T >
NumericVector<T>& NumericVector< T >::operator/= ( const T  a  )  [inline, inherited]

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

Definition at line 354 of file numeric_vector.h.

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

template<typename T >
NumericVector< T > & EpetraVector< T >::operator= ( const std::vector< T > &  v  )  [inline, 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 375 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, EpetraVector< T >::first_local_index(), EpetraVector< T >::local_size(), and EpetraVector< T >::size().

00376 {
00377   T * values = _vec->Values();
00378   
00383   if(this->size() == v.size())
00384   {
00385     const unsigned int nl=this->local_size();
00386     const unsigned int fli=this->first_local_index();
00387     
00388     for(unsigned int i=0;i<nl;i++)
00389       values[i]=v[fli+i];
00390   }
00391 
00396   else
00397   {
00398     libmesh_assert(v.size()==this->local_size());
00399     
00400     const unsigned int nl=this->local_size();
00401     
00402     for(unsigned int i=0;i<nl;i++)
00403       values[i]=v[i];
00404   }
00405 
00406   return *this;
00407 }

template<typename T >
EpetraVector< T > & EpetraVector< T >::operator= ( const EpetraVector< T > &  V  )  [inline]

$U = V$: copy all components.

Definition at line 364 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec.

00365 {
00366   (*_vec) = *v._vec;
00367   
00368   return *this;
00369 }

template<typename T >
NumericVector< T > & EpetraVector< T >::operator= ( const NumericVector< T > &  V  )  [inline, virtual]

$U = V$: copy all components.

Implements NumericVector< T >.

Definition at line 351 of file trilinos_epetra_vector.C.

00352 {
00353   const EpetraVector<T>* v = libmesh_cast_ptr<const EpetraVector<T>*>(&v_in);
00354 
00355   *this = *v;
00356   
00357   return *this;
00358 }

template<typename T >
NumericVector< T > & EpetraVector< T >::operator= ( const T  s  )  [inline, 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 340 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec.

00341 {
00342   _vec->PutScalar(s_in);
00343   
00344   return *this;
00345 }

template<typename T >
void EpetraVector< T >::pointwise_mult ( const NumericVector< T > &  vec1,
const NumericVector< T > &  vec2 
) [inline, 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 328 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec.

00330 {
00331   const EpetraVector<T>* V1 = libmesh_cast_ptr<const EpetraVector<T>*>(&vec1);
00332   const EpetraVector<T>* V2 = libmesh_cast_ptr<const EpetraVector<T>*>(&vec2);
00333 
00334   _vec->Multiply(1.0, *V1->_vec, *V2->_vec, 0.0);
00335 }

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

Definition at line 739 of file numeric_vector.h.

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

00740 {
00741   libmesh_assert (this->initialized());
00742   os << "Size\tglobal =  " << this->size()
00743      << "\t\tlocal =  " << this->local_size() << std::endl;
00744   
00745   // std::complex<>::operator<<() is defined, but use this form
00746   os << "#\tReal part\t\tImaginary part" << std::endl;
00747   for (unsigned int i=this->first_local_index(); i<this->last_local_index(); i++)
00748     os << i << "\t" 
00749        << (*this)(i).real() << "\t\t" 
00750        << (*this)(i).imag() << std::endl;
00751 }

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

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

00758 {
00759   libmesh_assert (this->initialized());
00760   os << "Size\tglobal =  " << this->size()
00761      << "\t\tlocal =  " << this->local_size() << std::endl;
00762 
00763   os << "#\tValue" << std::endl;
00764   for (unsigned int i=this->first_local_index(); i<this->last_local_index(); i++)
00765     os << i << "\t" << (*this)(i) << std::endl;
00766 }

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

Definition at line 772 of file numeric_vector.h.

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

00773 {
00774   libmesh_assert (this->initialized());
00775 
00776   std::vector<Complex> v(this->size());
00777   this->localize(v);
00778 
00779   // Right now we only want one copy of the output
00780   if (libMesh::processor_id())
00781     return;
00782   
00783   os << "Size\tglobal =  " << this->size() << std::endl;
00784   os << "#\tReal part\t\tImaginary part" << std::endl;
00785   for (unsigned int i=0; i!=v.size(); i++)
00786     os << i << "\t" 
00787        << v[i].real() << "\t\t" 
00788        << v[i].imag() << std::endl;
00789 }

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

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

00795 {
00796   libmesh_assert (this->initialized());
00797 
00798   std::vector<T> v(this->size());
00799   this->localize(v);
00800 
00801   // Right now we only want one copy of the output
00802   if (libMesh::processor_id())
00803     return;
00804 
00805   os << "Size\tglobal =  " << this->size() << std::endl;
00806   os << "#\tValue" << std::endl;
00807   for (unsigned int i=0; i!=v.size(); i++)
00808     os << i << "\t" << v[i] << std::endl;
00809 }

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().

00084 {
00085 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
00086   
00087   std::cout << ReferenceCounter::get_info();
00088   
00089 #endif
00090 }

template<typename T >
void EpetraVector< T >::print_matlab ( const std::string  name = "NULL"  )  const [inline, 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 584 of file trilinos_epetra_vector.C.

00585 {
00586   libmesh_not_implemented();
00587 
00588 //   libmesh_assert (this->initialized());
00589 //   libmesh_assert (this->closed());
00590   
00591 //   int ierr=0; 
00592 //   EpetraViewer epetra_viewer;
00593 
00594 
00595 //   ierr = EpetraViewerCreate (libMesh::COMM_WORLD,
00596 //                          &epetra_viewer);
00597 //          CHKERRABORT(libMesh::COMM_WORLD,ierr);
00598 
00599 //   /**
00600 //    * Create an ASCII file containing the matrix
00601 //    * if a filename was provided.  
00602 //    */
00603 //   if (name != "NULL")
00604 //     {
00605 //       ierr = EpetraViewerASCIIOpen( libMesh::COMM_WORLD,
00606 //                                 name.c_str(),
00607 //                                 &epetra_viewer);
00608 //              CHKERRABORT(libMesh::COMM_WORLD,ierr);
00609       
00610 //       ierr = EpetraViewerSetFormat (epetra_viewer,
00611 //                                 EPETRA_VIEWER_ASCII_MATLAB);
00612 //              CHKERRABORT(libMesh::COMM_WORLD,ierr);
00613   
00614 //       ierr = VecView (_vec, epetra_viewer);
00615 //              CHKERRABORT(libMesh::COMM_WORLD,ierr);
00616 //     }
00617 
00618 //   /**
00619 //    * Otherwise the matrix will be dumped to the screen.
00620 //    */
00621 //   else
00622 //     {
00623 //       ierr = EpetraViewerSetFormat (EPETRA_VIEWER_STDOUT_WORLD,
00624 //                                 EPETRA_VIEWER_ASCII_MATLAB);
00625 //              CHKERRABORT(libMesh::COMM_WORLD,ierr);
00626   
00627 //       ierr = VecView (_vec, EPETRA_VIEWER_STDOUT_WORLD);
00628 //              CHKERRABORT(libMesh::COMM_WORLD,ierr);
00629 //     }
00630 
00631 
00632 //   /**
00633 //    * Destroy the viewer.
00634 //    */
00635 //   ierr = EpetraViewerDestroy (epetra_viewer);
00636 //          CHKERRABORT(libMesh::COMM_WORLD,ierr);
00637 }

template<typename T >
int EpetraVector< T >::ReplaceGlobalValues ( int  numIDs,
const int *  GIDs,
const int *  numValuesPerID,
const double *  values 
) [inline, private]

Definition at line 795 of file trilinos_epetra_vector.C.

References EpetraVector< T >::inputValues().

00798 {
00799   return( inputValues( numIDs, GIDs, numValuesPerID, values, false) );
00800 }

template<typename T >
int EpetraVector< T >::ReplaceGlobalValues ( const Epetra_IntSerialDenseVector &  GIDs,
const Epetra_SerialDenseVector &  values 
) [inline, private]

Copy values into the vector, replacing any values that already exist for the specified GIDs.

Parameters:
GIDs List of global ids. Must be the same length as the accompanying list of values.
values List of coefficient values. Must be the same length as the accompanying list of GIDs.

Definition at line 783 of file trilinos_epetra_vector.C.

References EpetraVector< T >::inputValues().

00785 {
00786   if (GIDs.Length() != values.Length()) {
00787     return(-1);
00788   }
00789 
00790   return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), false) );
00791 }

template<typename T >
int EpetraVector< T >::ReplaceGlobalValues ( int  numIDs,
const int *  GIDs,
const double *  values 
) [inline, private]

Copy values into the vector overwriting any values that already exist for the specified indices.

Definition at line 775 of file trilinos_epetra_vector.C.

References EpetraVector< T >::inputValues().

Referenced by EpetraVector< T >::insert(), and EpetraVector< T >::set().

00777 {
00778   return( inputValues( numIDs, GIDs, values, false) );
00779 }

template<typename T >
void EpetraVector< T >::scale ( const T  factor  )  [inline, virtual]

Scale each element of the vector by the given factor.

Implements NumericVector< T >.

Definition at line 302 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec.

00303 {
00304   _vec->Scale(factor_in);
00305 }

template<typename T >
void EpetraVector< T >::set ( const unsigned int  i,
const T  value 
) [inline, virtual]

v(i) = value

Implements NumericVector< T >.

Definition at line 126 of file trilinos_epetra_vector.C.

References NumericVector< T >::_is_closed, EpetraVector< T >::ReplaceGlobalValues(), and EpetraVector< T >::size().

00127 {
00128   int i = static_cast<int> (i_in);
00129   T value = value_in;
00130 
00131   libmesh_assert(i_in<this->size());
00132 
00133   ReplaceGlobalValues(1, &i, &value);
00134 
00135   this->_is_closed = false;
00136 }

template<typename T >
void EpetraVector< T >::setIgnoreNonLocalEntries ( bool  flag  )  [inline, private]

Set whether or not non-local data values should be ignored.

Definition at line 557 of file trilinos_epetra_vector.h.

References EpetraVector< T >::ignoreNonLocalEntries_.

00557                                            {
00558     ignoreNonLocalEntries_ = flag;
00559   }

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

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

Implements NumericVector< T >.

Definition at line 865 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec, and NumericVector< T >::initialized().

Referenced by EpetraVector< T >::add(), EpetraVector< T >::localize(), EpetraVector< T >::localize_to_one(), EpetraVector< T >::operator=(), and EpetraVector< T >::set().

00866 {
00867   libmesh_assert (this->initialized());  
00868 
00869   return _vec->GlobalLength();
00870 }

template<class T >
Real NumericVector< T >::subset_l1_norm ( const std::set< unsigned int > &  indices  )  [inline, 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.

References NumericVector< T >::sum().

Referenced by System::discrete_var_norm().

00194 {
00195   NumericVector<T> & v = *this;
00196   
00197   std::set<unsigned int>::iterator it = indices.begin();
00198   const std::set<unsigned int>::iterator it_end = indices.end();
00199 
00200   Real norm = 0;
00201   
00202   for(; it!=it_end; ++it)
00203     norm += std::abs(v(*it));
00204 
00205   Parallel::sum(norm);
00206 
00207   return norm;
00208 }

template<class T >
Real NumericVector< T >::subset_l2_norm ( const std::set< unsigned int > &  indices  )  [inline, 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(), and NumericVector< T >::sum().

Referenced by System::discrete_var_norm().

00212 {
00213   NumericVector<T> & v = *this;
00214 
00215   std::set<unsigned int>::iterator it = indices.begin();
00216   const std::set<unsigned int>::iterator it_end = indices.end();
00217 
00218   Real norm = 0;
00219   
00220   for(; it!=it_end; ++it)
00221     norm += libmesh_norm(v(*it));
00222 
00223   Parallel::sum(norm);
00224 
00225   return std::sqrt(norm);
00226 }

template<class T >
Real NumericVector< T >::subset_linfty_norm ( const std::set< unsigned int > &  indices  )  [inline, 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 NumericVector< T >::abs(), and NumericVector< T >::max().

Referenced by System::discrete_var_norm().

00230 {
00231   NumericVector<T> & v = *this;
00232 
00233   std::set<unsigned int>::iterator it = indices.begin();
00234   const std::set<unsigned int>::iterator it_end = indices.end();
00235 
00236   Real norm = 0;
00237   
00238   for(; it!=it_end; ++it)
00239     {
00240       Real value = std::abs(v(*it));
00241       if(value > norm)
00242         norm = value;
00243     }
00244 
00245   Parallel::max(norm);
00246 
00247   return norm;
00248 }

template<typename T >
T EpetraVector< T >::sum (  )  const [inline, virtual]

Returns:
the sum of values in a vector

Implements NumericVector< T >.

Definition at line 45 of file trilinos_epetra_vector.C.

References EpetraVector< T >::_vec, and NumericVector< T >::closed().

00046 {
00047   libmesh_assert(this->closed());
00048   
00049   const unsigned int nl = _vec->MyLength();
00050 
00051   T sum=0.0;
00052 
00053   T * values = _vec->Values();
00054   
00055   for (unsigned int i=0; i<nl; i++)
00056     sum += values[i];
00057 
00058   Parallel::sum<T>(sum);
00059   
00060   return sum;
00061 }

template<typename T >
int EpetraVector< T >::SumIntoGlobalValues ( int  numIDs,
const int *  GIDs,
const int *  numValuesPerID,
const double *  values 
) [inline, private]

Definition at line 766 of file trilinos_epetra_vector.C.

References EpetraVector< T >::inputValues().

00769 {
00770   return( inputValues( numIDs, GIDs, numValuesPerID, values, true) );
00771 }

template<typename T >
int EpetraVector< T >::SumIntoGlobalValues ( const Epetra_IntSerialDenseVector &  GIDs,
const Epetra_SerialDenseVector &  values 
) [inline, private]

Accumulate values into the vector, adding them to any values that already exist for the specified GIDs.

Parameters:
GIDs List of global ids. Must be the same length as the accompanying list of values.
values List of coefficient values. Must be the same length as the accompanying list of GIDs.

Definition at line 754 of file trilinos_epetra_vector.C.

References EpetraVector< T >::inputValues().

00756 {
00757   if (GIDs.Length() != values.Length()) {
00758     return(-1);
00759   }
00760 
00761   return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), true) );
00762 }

template<typename T >
int EpetraVector< T >::SumIntoGlobalValues ( int  numIDs,
const int *  GIDs,
const double *  values 
) [inline, private]

Accumulate values into the vector, adding them to any values that already exist for the specified indices.

Definition at line 746 of file trilinos_epetra_vector.C.

References EpetraVector< T >::inputValues().

Referenced by EpetraVector< T >::add(), and EpetraVector< T >::add_vector().

00748 {
00749   return( inputValues( numIDs, GIDs, values, true) );
00750 }

template<typename T >
void EpetraVector< T >::swap ( NumericVector< T > &  v  )  [inline, virtual]

Swaps the raw Epetra vector context pointers.

Reimplemented from NumericVector< T >.

Definition at line 949 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_destroy_vec_on_exit, and EpetraVector< T >::_vec.

00950 {
00951   EpetraVector<T>& v = libmesh_cast_ref<EpetraVector<T>&>(other);
00952 
00953   std::swap(_vec, v._vec);
00954   std::swap(_destroy_vec_on_exit, v._destroy_vec_on_exit);
00955 }

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.

00119 { return _type; }

template<typename T >
Epetra_Vector* EpetraVector< 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 474 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec.

Referenced by EpetraMatrix< T >::get_diagonal(), and NoxNonlinearSolver< T >::solve().

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

template<typename T >
void EpetraVector< T >::zero (  )  [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 826 of file trilinos_epetra_vector.h.

References EpetraVector< T >::_vec, and NumericVector< T >::initialized().

Referenced by Problem_Interface::computeF(), and EpetraVector< T >::init().

00827 {
00828   libmesh_assert (this->initialized());
00829 
00830   _vec->PutScalar(0.0);
00831 }

template<typename T >
AutoPtr< NumericVector< T > > EpetraVector< T >::zero_clone (  )  const [inline, virtual]

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

Implements NumericVector< T >.

Definition at line 837 of file trilinos_epetra_vector.h.

00838 {
00839   AutoPtr<NumericVector<T> > cloned_vector (new EpetraVector<T>);
00840 
00841   cloned_vector->init(*this);
00842 
00843   return cloned_vector;
00844 }


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

00546   {
00547     v.print_global(os);
00548     return os;
00549   }


Member Data Documentation

template<typename T >
bool EpetraVector< T >::_destroy_vec_on_exit [private]

This boolean value should only be set to false for the constructor which takes a Epetra Vec object.

Definition at line 493 of file trilinos_epetra_vector.h.

Referenced by EpetraVector< T >::clear(), and EpetraVector< T >::swap().

template<typename T >
Epetra_Map* EpetraVector< T >::_map [private]

Holds the distributed Map

Definition at line 487 of file trilinos_epetra_vector.h.

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

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 >
double* EpetraVector< T >::myCoefs_ [private]

template<typename T >
int EpetraVector< T >::myFirstID_ [private]

template<typename T >
int EpetraVector< T >::myNumIDs_ [private]


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

Site Created By: libMesh Developers
Last modified: November 25 2009 03:44:10.

Hosted By:
SourceForge.net Logo