libMesh::EpetraVector< T > Class Template Reference
#include <trilinos_epetra_vector.h>

Public Member Functions | |
| EpetraVector (const ParallelType type=AUTOMATIC) | |
| EpetraVector (const numeric_index_type n, const ParallelType type=AUTOMATIC) | |
| EpetraVector (const numeric_index_type n, const numeric_index_type n_local, const ParallelType type=AUTOMATIC) | |
| EpetraVector (const numeric_index_type N, const numeric_index_type n_local, const std::vector< numeric_index_type > &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 numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType type=AUTOMATIC) |
| void | init (const numeric_index_type N, const bool fast=false, const ParallelType type=AUTOMATIC) |
| void | init (const numeric_index_type, const numeric_index_type, const std::vector< numeric_index_type > &, 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 |
| T | sum () const |
| Real | l1_norm () const |
| Real | l2_norm () const |
| Real | linfty_norm () const |
| numeric_index_type | size () const |
| numeric_index_type | local_size () const |
| numeric_index_type | first_local_index () const |
| numeric_index_type | last_local_index () const |
| T | operator() (const numeric_index_type i) const |
| NumericVector< T > & | operator+= (const NumericVector< T > &V) |
| NumericVector< T > & | operator-= (const NumericVector< T > &V) |
| virtual void | reciprocal () |
| void | set (const numeric_index_type i, const T value) |
| void | add (const numeric_index_type 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< numeric_index_type > &dof_indices) |
| void | add_vector (const NumericVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| void | add_vector (const NumericVector< T > &V, const SparseMatrix< T > &A) |
| void | add_vector (const DenseVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| void | add_vector_transpose (const NumericVector< T > &V, const SparseMatrix< T > &A_trans) |
| virtual void | insert (const std::vector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
| virtual void | insert (const NumericVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| virtual void | insert (const DenseVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| virtual void | insert (const DenseSubVector< T > &V, const std::vector< numeric_index_type > &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< numeric_index_type > &send_list) const |
| void | localize (const numeric_index_type first_local_idx, const numeric_index_type last_local_idx, const std::vector< numeric_index_type > &send_list) |
| void | localize_to_one (std::vector< T > &v_local, const processor_id_type 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< numeric_index_type > &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 void | init (const NumericVector< T > &other, const bool fast=false)=0 |
| virtual Real | subset_l1_norm (const std::set< numeric_index_type > &indices) const |
| virtual Real | subset_l2_norm (const std::set< numeric_index_type > &indices) const |
| virtual Real | subset_linfty_norm (const std::set< numeric_index_type > &indices) const |
| virtual T | el (const numeric_index_type i) const |
| virtual void | get (const std::vector< numeric_index_type > &index, std::vector< T > &values) const |
| virtual NumericVector< T > & | operator+= (const NumericVector< T > &V)=0 |
| virtual NumericVector< T > & | operator-= (const NumericVector< T > &V)=0 |
| NumericVector< T > & | operator*= (const T a) |
| NumericVector< T > & | operator/= (const T a) |
| virtual void | add (const NumericVector< T > &V)=0 |
| virtual void | add (const T a, const NumericVector< T > &v)=0 |
| virtual void | add_vector (const NumericVector< T > &V, const std::vector< numeric_index_type > &dof_indices)=0 |
| virtual void | add_vector (const NumericVector< T > &, const SparseMatrix< T > &)=0 |
| void | add_vector (const NumericVector< T > &v, const ShellMatrix< T > &a) |
| virtual void | add_vector_transpose (const NumericVector< T > &, const SparseMatrix< T > &)=0 |
| virtual void | insert (const NumericVector< T > &V, const std::vector< numeric_index_type > &dof_indices)=0 |
| virtual T | dot (const NumericVector< T > &) const =0 |
| virtual void | localize (NumericVector< T > &v_local) const =0 |
| virtual void | localize (NumericVector< T > &v_local, const std::vector< numeric_index_type > &send_list) const =0 |
| virtual int | compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
| virtual int | local_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
| virtual int | global_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
| virtual void | pointwise_mult (const NumericVector< T > &vec1, const NumericVector< T > &vec2)=0 |
| virtual void | print (std::ostream &os=libMesh::out) const |
| template<> | |
| void | print (std::ostream &os) const |
| virtual void | print_global (std::ostream &os=libMesh::out) const |
| template<> | |
| void | print_global (std::ostream &os) const |
| virtual void | create_subvector (NumericVector< T > &, const std::vector< numeric_index_type > &) const |
| virtual void | swap (NumericVector< T > &v) |
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 (std::ostream &out=libMesh::out) |
| static unsigned int | n_objects () |
| static void | enable_print_counter_info () |
| static void | disable_print_counter_info () |
Protected Types | |
| typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
| void | increment_constructor_count (const std::string &name) |
| void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
| bool | _is_closed |
| bool | _is_initialized |
| ParallelType | _type |
Static Protected Attributes | |
| static Counts | _counts |
| static Threads::atomic < unsigned int > | _n_objects |
| static Threads::spin_mutex | _mutex |
| static bool | _enable_print_counter = true |
Private 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 libMesh::EpetraVector< T >
Epetra vector. Provides a nice interface to the Trilinos Epetra data structures for parallel vectors.
Definition at line 59 of file trilinos_epetra_vector.h.
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 113 of file reference_counter.h.
Constructor & Destructor Documentation
| libMesh::EpetraVector< T >::EpetraVector | ( | const ParallelType | type = AUTOMATIC |
) | [inline, explicit] |
Dummy-Constructor. Dimension=0
Definition at line 615 of file trilinos_epetra_vector.h.
References libMesh::NumericVector< T >::_type.
00616 : _destroy_vec_on_exit(true), 00617 myFirstID_(0), 00618 myNumIDs_(0), 00619 myCoefs_(NULL), 00620 nonlocalIDs_(NULL), 00621 nonlocalElementSize_(NULL), 00622 numNonlocalIDs_(0), 00623 allocatedNonlocalLength_(0), 00624 nonlocalCoefs_(NULL), 00625 ignoreNonLocalEntries_(false) 00626 { 00627 this->_type = type; 00628 }
| libMesh::EpetraVector< T >::EpetraVector | ( | const numeric_index_type | n, | |
| const ParallelType | type = AUTOMATIC | |||
| ) | [inline, explicit] |
Constructor. Set dimension to n and initialize all elements with zero.
Definition at line 634 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
00636 : _destroy_vec_on_exit(true), 00637 myFirstID_(0), 00638 myNumIDs_(0), 00639 myCoefs_(NULL), 00640 nonlocalIDs_(NULL), 00641 nonlocalElementSize_(NULL), 00642 numNonlocalIDs_(0), 00643 allocatedNonlocalLength_(0), 00644 nonlocalCoefs_(NULL), 00645 ignoreNonLocalEntries_(false) 00646 00647 { 00648 this->init(n, n, false, type); 00649 }
| libMesh::EpetraVector< T >::EpetraVector | ( | const numeric_index_type | n, | |
| const numeric_index_type | 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 655 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
00658 : _destroy_vec_on_exit(true), 00659 myFirstID_(0), 00660 myNumIDs_(0), 00661 myCoefs_(NULL), 00662 nonlocalIDs_(NULL), 00663 nonlocalElementSize_(NULL), 00664 numNonlocalIDs_(0), 00665 allocatedNonlocalLength_(0), 00666 nonlocalCoefs_(NULL), 00667 ignoreNonLocalEntries_(false) 00668 { 00669 this->init(n, n_local, false, type); 00670 }
| libMesh::EpetraVector< T >::EpetraVector | ( | const numeric_index_type | N, | |
| const numeric_index_type | n_local, | |||
| const std::vector< numeric_index_type > & | 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 709 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
00713 : _destroy_vec_on_exit(true), 00714 myFirstID_(0), 00715 myNumIDs_(0), 00716 myCoefs_(NULL), 00717 nonlocalIDs_(NULL), 00718 nonlocalElementSize_(NULL), 00719 numNonlocalIDs_(0), 00720 allocatedNonlocalLength_(0), 00721 nonlocalCoefs_(NULL), 00722 ignoreNonLocalEntries_(false) 00723 { 00724 this->init(n, n_local, ghost, false, type); 00725 }
| libMesh::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 677 of file trilinos_epetra_vector.h.
References libMesh::NumericVector< T >::_is_closed, libMesh::NumericVector< T >::_is_initialized, libMesh::NumericVector< T >::_type, libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::myCoefs_, libMesh::EpetraVector< T >::myFirstID_, libMesh::EpetraVector< T >::myNumIDs_, and libMeshEnums::PARALLEL.
00678 : _destroy_vec_on_exit(false), 00679 myFirstID_(0), 00680 myNumIDs_(0), 00681 myCoefs_(NULL), 00682 nonlocalIDs_(NULL), 00683 nonlocalElementSize_(NULL), 00684 numNonlocalIDs_(0), 00685 allocatedNonlocalLength_(0), 00686 nonlocalCoefs_(NULL), 00687 ignoreNonLocalEntries_(false) 00688 { 00689 _vec = &v; 00690 00691 this->_type = PARALLEL; // FIXME - need to determine this from v! 00692 00693 myFirstID_ = _vec->Map().MinMyGID(); 00694 myNumIDs_ = _vec->Map().NumMyElements(); 00695 00696 //Currently we impose the restriction that NumVectors==1, so we won't 00697 //need the LDA argument when calling ExtractView. Hence the "dummy" arg. 00698 int dummy; 00699 _vec->ExtractView(&myCoefs_, &dummy); 00700 00701 this->_is_closed = true; 00702 this->_is_initialized = true; 00703 }
| libMesh::EpetraVector< T >::~EpetraVector | ( | ) | [inline] |
Destructor, deallocates memory. Made virtual to allow for derived classes to behave properly.
Definition at line 742 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::clear().
00743 { 00744 this->clear (); 00745 }
Member Function Documentation
| void libMesh::EpetraVector< T >::abs | ( | ) | [inline, virtual] |
v = abs(v)... that is, each entry in v is replaced by its absolute value.
Implements libMesh::NumericVector< T >.
Definition at line 346 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
| virtual void libMesh::NumericVector< T >::add | ( | const T | a, | |
| const NumericVector< T > & | v | |||
| ) | [pure virtual, inherited] |
. Simple vector addition, equal to the operator +=.
| virtual void libMesh::NumericVector< T >::add | ( | const NumericVector< T > & | V | ) | [pure virtual, inherited] |
: Simple vector addition, equal to the operator +=.
| void libMesh::EpetraVector< T >::add | ( | const T | a, | |
| const NumericVector< T > & | v | |||
| ) | [inline] |
. Simple vector addition, equal to the operator +=.
Definition at line 275 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, and libMesh::EpetraVector< T >::size().
| void libMesh::EpetraVector< T >::add | ( | const NumericVector< T > & | V | ) | [inline] |
. Simple vector addition, equal to the operator +=.
Definition at line 268 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::add().
00269 { 00270 this->add (1., v); 00271 }
| void libMesh::EpetraVector< T >::add | ( | const T | s | ) | [inline, virtual] |
. Addition of s to all components. Note that s is a scalar and not a vector.
Implements libMesh::NumericVector< T >.
Definition at line 254 of file trilinos_epetra_vector.C.
References libMesh::NumericVector< T >::_is_closed, and libMesh::EpetraVector< T >::_vec.
00255 { 00256 const unsigned int nl = _vec->MyLength(); 00257 00258 T * values = _vec->Values(); 00259 00260 for (unsigned int i=0; i<nl; i++) 00261 values[i]+=v_in; 00262 00263 this->_is_closed = false; 00264 }
| void libMesh::EpetraVector< T >::add | ( | const numeric_index_type | i, | |
| const T | value | |||
| ) | [virtual] |
v(i) += value
Implements libMesh::NumericVector< T >.
Referenced by libMesh::EpetraVector< T >::add(), libMesh::EpetraVector< T >::operator+=(), and libMesh::EpetraVector< T >::operator-=().
| void libMesh::NumericVector< T >::add_vector | ( | const NumericVector< T > & | v, | |
| const ShellMatrix< T > & | a | |||
| ) | [inline, inherited] |
, add the product of a ShellMatrix A and a NumericVector V to this, where this=U.
Definition at line 358 of file numeric_vector.C.
References libMesh::ShellMatrix< T >::vector_mult_add().
| virtual void libMesh::NumericVector< T >::add_vector | ( | const NumericVector< T > & | , | |
| const SparseMatrix< T > & | ||||
| ) | [pure virtual, inherited] |
, add the product of a SparseMatrix A and a NumericVector V to this, where this=U.
| virtual void libMesh::NumericVector< T >::add_vector | ( | const NumericVector< T > & | V, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [pure virtual, inherited] |
, where U and V are type NumericVector<T> and you want to specify WHERE to add the NumericVector<T> V
| void libMesh::EpetraVector< T >::add_vector | ( | const DenseVector< T > & | V, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [virtual] |
where U and V are type DenseVector<T> and you want to specify WHERE to add the DenseVector<T> V
Implements libMesh::NumericVector< T >.
| void libMesh::EpetraVector< T >::add_vector | ( | const NumericVector< T > & | V, | |
| const SparseMatrix< T > & | A | |||
| ) | [inline] |
, add the product of a SparseMatrix A and a NumericVector V to this, where this=U.
Definition at line 216 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, and libMesh::EpetraVector< T >::zero_clone().
00218 { 00219 const EpetraVector<T>* V = libmesh_cast_ptr<const EpetraVector<T>*>(&V_in); 00220 const EpetraMatrix<T>* A = libmesh_cast_ptr<const EpetraMatrix<T>*>(&A_in); 00221 00222 // FIXME - does Trilinos let us do this *without* memory allocation? 00223 AutoPtr<NumericVector<T> > temp = V->zero_clone(); 00224 EpetraVector<T>* tempV = libmesh_cast_ptr<EpetraVector<T>*>(temp.get()); 00225 A->mat()->Multiply(false, *V->_vec, *tempV->_vec); 00226 *this += *temp; 00227 }
| void libMesh::EpetraVector< T >::add_vector | ( | const NumericVector< T > & | V, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) |
where U and V are type NumericVector<T> and you want to specify WHERE to add the NumericVector<T> V
| void libMesh::EpetraVector< T >::add_vector | ( | const std::vector< T > & | v, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [virtual] |
where v is a std::vector<T> and you want to specify WHERE to add it
Implements libMesh::NumericVector< T >.
| virtual void libMesh::NumericVector< T >::add_vector_transpose | ( | const NumericVector< T > & | , | |
| const SparseMatrix< T > & | ||||
| ) | [pure virtual, inherited] |
, add the product of the transpose of a SparseMatrix A_trans and a NumericVector V to this, where this=U.
| void libMesh::EpetraVector< T >::add_vector_transpose | ( | const NumericVector< T > & | V, | |
| const SparseMatrix< T > & | A_trans | |||
| ) | [inline] |
, add the product of the transpose of a SparseMatrix A_trans and a NumericVector V to this, where this=U.
Definition at line 245 of file trilinos_epetra_vector.C.
| AutoPtr< NumericVector< T > > libMesh::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 45 of file numeric_vector.C.
References libMesh::LASPACK_SOLVERS, libMeshEnums::PETSC_SOLVERS, and libMesh::TRILINOS_SOLVERS.
Referenced by libMesh::ExactErrorEstimator::estimate_error().
00046 { 00047 // Build the appropriate vector 00048 switch (solver_package) 00049 { 00050 00051 00052 #ifdef LIBMESH_HAVE_LASPACK 00053 case LASPACK_SOLVERS: 00054 { 00055 AutoPtr<NumericVector<T> > ap(new LaspackVector<T>); 00056 return ap; 00057 } 00058 #endif 00059 00060 00061 #ifdef LIBMESH_HAVE_PETSC 00062 case PETSC_SOLVERS: 00063 { 00064 AutoPtr<NumericVector<T> > ap(new PetscVector<T>); 00065 return ap; 00066 } 00067 #endif 00068 00069 00070 #ifdef LIBMESH_HAVE_TRILINOS 00071 case TRILINOS_SOLVERS: 00072 { 00073 AutoPtr<NumericVector<T> > ap(new EpetraVector<T>); 00074 return ap; 00075 } 00076 #endif 00077 00078 00079 default: 00080 AutoPtr<NumericVector<T> > ap(new DistributedVector<T>); 00081 return ap; 00082 } 00083 00084 AutoPtr<NumericVector<T> > ap(NULL); 00085 return ap; 00086 }
| void libMesh::EpetraVector< T >::clear | ( | ) | [inline, virtual] |
- Returns:
- the
EpetraVector<T>to a pristine state.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 843 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_destroy_vec_on_exit, libMesh::NumericVector< T >::_is_closed, libMesh::NumericVector< T >::_is_initialized, libMesh::EpetraVector< T >::_vec, and libMesh::NumericVector< T >::initialized().
Referenced by libMesh::EpetraVector< T >::~EpetraVector().
00844 { 00845 if ((this->initialized()) && (this->_destroy_vec_on_exit)) 00846 { 00847 delete _vec; 00848 _vec = NULL; 00849 } 00850 00851 this->_is_closed = this->_is_initialized = false; 00852 }
| AutoPtr< NumericVector< T > > libMesh::EpetraVector< T >::clone | ( | ) | const [inline, virtual] |
Creates a copy of this vector and returns it in an AutoPtr.
Implements libMesh::NumericVector< T >.
Definition at line 883 of file trilinos_epetra_vector.h.
| void libMesh::EpetraVector< T >::close | ( | ) | [inline, virtual] |
Call the assemble functions
Implements libMesh::NumericVector< T >.
Definition at line 830 of file trilinos_epetra_vector.h.
References libMesh::NumericVector< T >::_is_closed, libMesh::EpetraVector< T >::GlobalAssemble(), and libMesh::NumericVector< T >::initialized().
Referenced by libMesh::Problem_Interface::computeF(), and libMesh::EpetraVector< T >::reciprocal().
00831 { 00832 libmesh_assert (this->initialized()); 00833 00834 this->GlobalAssemble(); 00835 00836 this->_is_closed = true; 00837 }
| virtual bool libMesh::NumericVector< T >::closed | ( | ) | const [inline, virtual, inherited] |
- Returns:
- true if the vector is closed and ready for computation, false otherwise.
Definition at line 127 of file numeric_vector.h.
Referenced by libMesh::DofMap::enforce_constraints_exactly(), libMesh::EpetraVector< T >::l1_norm(), libMesh::PetscVector< T >::l1_norm(), libMesh::LaspackVector< T >::l1_norm(), libMesh::EpetraVector< T >::l2_norm(), libMesh::PetscVector< T >::l2_norm(), libMesh::LaspackVector< T >::l2_norm(), libMesh::EpetraVector< T >::linfty_norm(), libMesh::PetscVector< T >::linfty_norm(), libMesh::LaspackVector< T >::linfty_norm(), libMesh::DofMap::max_constraint_error(), libMesh::EpetraVector< T >::operator+=(), libMesh::PetscVector< T >::operator+=(), libMesh::LaspackVector< T >::operator+=(), libMesh::DistributedVector< T >::operator+=(), libMesh::EpetraVector< T >::operator-=(), libMesh::PetscVector< T >::operator-=(), libMesh::LaspackVector< T >::operator-=(), libMesh::DistributedVector< T >::operator-=(), libMesh::PetscVector< T >::operator=(), libMesh::LaspackVector< T >::operator=(), libMesh::PetscVector< T >::print_matlab(), libMesh::EpetraVector< T >::sum(), libMesh::PetscVector< T >::sum(), libMesh::LaspackVector< T >::sum(), libMesh::EpetraVector< T >::zero(), libMesh::PetscVector< T >::zero(), and libMesh::LaspackVector< T >::zero().
00127 { return _is_closed; }
| int libMesh::NumericVector< T >::compare | ( | const NumericVector< T > & | other_vector, | |
| const Real | threshold = TOLERANCE | |||
| ) | const [inline, virtual, inherited] |
- Returns:
-1whenthisis equivalent toother_vector, up to the giventhreshold. When differences occur, the return value contains the first indexiwhere the difference(a[i]-b[i]) exceeded the threshold. When no threshold is given, thelibMeshTOLERANCEis used.
Definition at line 90 of file numeric_vector.C.
References std::abs(), libMesh::CommWorld, libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::last_local_index(), std::max(), libMesh::NumericVector< T >::max(), and libMesh::Parallel::Communicator::min().
00092 { 00093 libmesh_assert (this->initialized()); 00094 libmesh_assert (other_vector.initialized()); 00095 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00096 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00097 00098 int first_different_i = std::numeric_limits<int>::max(); 00099 numeric_index_type i = first_local_index(); 00100 00101 do 00102 { 00103 if ( std::abs( (*this)(i) - other_vector(i) ) > threshold ) 00104 first_different_i = i; 00105 else 00106 i++; 00107 } 00108 while (first_different_i==std::numeric_limits<int>::max() 00109 && i<last_local_index()); 00110 00111 // Find the correct first differing index in parallel 00112 CommWorld.min(first_different_i); 00113 00114 if (first_different_i == std::numeric_limits<int>::max()) 00115 return -1; 00116 00117 return first_different_i; 00118 }
| virtual void libMesh::NumericVector< T >::create_subvector | ( | NumericVector< T > & | , | |
| const std::vector< numeric_index_type > & | ||||
| ) | const [inline, virtual, inherited] |
Creates the subvector "subvector" from the indices in the "rows" array. Similar to the create_submatrix routine for the SparseMatrix class, it is currently only implemented for PetscVectors.
Definition at line 608 of file numeric_vector.h.
00610 { 00611 libMesh::err << "ERROR: Not Implemented in base class yet!" << std::endl; 00612 libmesh_error(); 00613 }
| virtual void libMesh::EpetraVector< T >::create_subvector | ( | NumericVector< T > & | subvector, | |
| const std::vector< numeric_index_type > & | rows | |||
| ) | const [virtual] |
Creates a "subvector" from this vector using the rows indices of the "rows" array.
| void libMesh::EpetraVector< T >::destroyNonlocalData | ( | ) | [inline, private] |
Definition at line 1078 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::allocatedNonlocalLength_, libMesh::EpetraVector< T >::nonlocalCoefs_, libMesh::EpetraVector< T >::nonlocalElementSize_, libMesh::EpetraVector< T >::nonlocalIDs_, and libMesh::EpetraVector< T >::numNonlocalIDs_.
Referenced by libMesh::EpetraVector< T >::FEoperatorequals(), and libMesh::EpetraVector< T >::GlobalAssemble().
01079 { 01080 if (allocatedNonlocalLength_ > 0) { 01081 delete [] nonlocalIDs_; 01082 delete [] nonlocalElementSize_; 01083 nonlocalIDs_ = NULL; 01084 nonlocalElementSize_ = NULL; 01085 for(int i=0; i<numNonlocalIDs_; ++i) { 01086 delete [] nonlocalCoefs_[i]; 01087 } 01088 delete [] nonlocalCoefs_; 01089 nonlocalCoefs_ = NULL; 01090 numNonlocalIDs_ = 0; 01091 allocatedNonlocalLength_ = 0; 01092 } 01093 return; 01094 }
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00107 { 00108 _enable_print_counter = false; 00109 return; 00110 }
| virtual T libMesh::NumericVector< T >::dot | ( | const NumericVector< T > & | ) | const [pure virtual, inherited] |
Computes the dot product, p = U.V
Referenced by libMesh::ContinuationSystem::continuation_solve(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ContinuationSystem::solve_tangent(), and libMesh::ContinuationSystem::update_solution().
| T libMesh::EpetraVector< T >::dot | ( | const NumericVector< T > & | V | ) | const [inline, virtual] |
Computes the dot product, p = U.V
Definition at line 353 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
00354 { 00355 const EpetraVector<T>* V = libmesh_cast_ptr<const EpetraVector<T>*>(&V_in); 00356 00357 T result=0.0; 00358 00359 _vec->Dot(*V->_vec, &result); 00360 00361 return result; 00362 }
| virtual T libMesh::NumericVector< T >::el | ( | const numeric_index_type | i | ) | const [inline, virtual, inherited] |
- Returns:
- the element
U(i)
Definition at line 324 of file numeric_vector.h.
| void libMesh::ReferenceCounter::enable_print_counter_info | ( | ) | [static, inherited] |
Methods to enable/disable the reference counter output from print_info()
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00101 { 00102 _enable_print_counter = true; 00103 return; 00104 }
| void libMesh::EpetraVector< T >::FEoperatorequals | ( | const EpetraVector< T > & | source | ) | [inline, private] |
Definition at line 1051 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::allocatedNonlocalLength_, libMesh::EpetraVector< T >::destroyNonlocalData(), libMesh::EpetraVector< T >::nonlocalCoefs_, libMesh::EpetraVector< T >::nonlocalElementSize_, libMesh::EpetraVector< T >::nonlocalIDs_, and libMesh::EpetraVector< T >::numNonlocalIDs_.
01052 { 01053 (*_vec) = *(source._vec); 01054 01055 destroyNonlocalData(); 01056 01057 if (source.allocatedNonlocalLength_ > 0) { 01058 allocatedNonlocalLength_ = source.allocatedNonlocalLength_; 01059 numNonlocalIDs_ = source.numNonlocalIDs_; 01060 nonlocalIDs_ = new int[allocatedNonlocalLength_]; 01061 nonlocalElementSize_ = new int[allocatedNonlocalLength_]; 01062 nonlocalCoefs_ = new double*[allocatedNonlocalLength_]; 01063 for(int i=0; i<numNonlocalIDs_; ++i) { 01064 int elemSize = source.nonlocalElementSize_[i]; 01065 nonlocalCoefs_[i] = new double[elemSize]; 01066 nonlocalIDs_[i] = source.nonlocalIDs_[i]; 01067 nonlocalElementSize_[i] = elemSize; 01068 for(int j=0; j<elemSize; ++j) { 01069 nonlocalCoefs_[i][j] = source.nonlocalCoefs_[i][j]; 01070 } 01071 } 01072 } 01073 }
| numeric_index_type libMesh::EpetraVector< T >::first_local_index | ( | ) | const [inline, virtual] |
- Returns:
- the index of the first vector element actually stored on this processor
Implements libMesh::NumericVector< T >.
Definition at line 918 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, and libMesh::NumericVector< T >::initialized().
Referenced by libMesh::EpetraVector< T >::operator()(), and libMesh::EpetraVector< T >::operator=().
00919 { 00920 libmesh_assert (this->initialized()); 00921 00922 return _vec->Map().MinMyGID(); 00923 }
| void libMesh::NumericVector< T >::get | ( | const std::vector< numeric_index_type > & | 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 libMesh::PetscVector< T >.
Definition at line 760 of file numeric_vector.h.
| std::string libMesh::ReferenceCounter::get_info | ( | ) | [static, inherited] |
Gets a string containing the reference information.
Definition at line 47 of file reference_counter.C.
References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().
Referenced by libMesh::ReferenceCounter::print_info().
00048 { 00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00050 00051 std::ostringstream oss; 00052 00053 oss << '\n' 00054 << " ---------------------------------------------------------------------------- \n" 00055 << "| Reference count information |\n" 00056 << " ---------------------------------------------------------------------------- \n"; 00057 00058 for (Counts::iterator it = _counts.begin(); 00059 it != _counts.end(); ++it) 00060 { 00061 const std::string name(it->first); 00062 const unsigned int creations = it->second.first; 00063 const unsigned int destructions = it->second.second; 00064 00065 oss << "| " << name << " reference count information:\n" 00066 << "| Creations: " << creations << '\n' 00067 << "| Destructions: " << destructions << '\n'; 00068 } 00069 00070 oss << " ---------------------------------------------------------------------------- \n"; 00071 00072 return oss.str(); 00073 00074 #else 00075 00076 return ""; 00077 00078 #endif 00079 }
| int libMesh::NumericVector< T >::global_relative_compare | ( | const NumericVector< T > & | other_vector, | |
| const Real | threshold = TOLERANCE | |||
| ) | const [inline, virtual, inherited] |
- Returns:
-1whenthisis equivalent toother_vector, up to the given local relativethreshold. When differences occur, the return value contains the first index where the difference(a[i]-b[i])/max_j(a[j],b[j]) exceeded the threshold. When no threshold is given, thelibMeshTOLERANCEis used.
Definition at line 155 of file numeric_vector.C.
References std::abs(), libMesh::CommWorld, libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::NumericVector< T >::linfty_norm(), std::max(), libMesh::NumericVector< T >::max(), libMesh::Parallel::Communicator::min(), and libMesh::Real.
00157 { 00158 libmesh_assert (this->initialized()); 00159 libmesh_assert (other_vector.initialized()); 00160 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00161 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00162 00163 int first_different_i = std::numeric_limits<int>::max(); 00164 numeric_index_type i = first_local_index(); 00165 00166 const Real my_norm = this->linfty_norm(); 00167 const Real other_norm = other_vector.linfty_norm(); 00168 const Real abs_threshold = std::max(my_norm, other_norm) * threshold; 00169 00170 do 00171 { 00172 if ( std::abs( (*this)(i) - other_vector(i) ) > abs_threshold ) 00173 first_different_i = i; 00174 else 00175 i++; 00176 } 00177 while (first_different_i==std::numeric_limits<int>::max() 00178 && i<last_local_index()); 00179 00180 // Find the correct first differing index in parallel 00181 CommWorld.min(first_different_i); 00182 00183 if (first_different_i == std::numeric_limits<int>::max()) 00184 return -1; 00185 00186 return first_different_i; 00187 }
| int libMesh::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 1006 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::destroyNonlocalData(), libMesh::EpetraVector< T >::ignoreNonLocalEntries_, libMesh::EpetraVector< T >::nonlocalCoefs_, libMesh::EpetraVector< T >::nonlocalElementSize_, libMesh::EpetraVector< T >::nonlocalIDs_, and libMesh::EpetraVector< T >::numNonlocalIDs_.
Referenced by libMesh::EpetraVector< T >::close().
01007 { 01008 //In this method we need to gather all the non-local (overlapping) data 01009 //that's been input on each processor, into the (probably) non-overlapping 01010 //distribution defined by the map that 'this' vector was constructed with. 01011 01012 //We don't need to do anything if there's only one processor or if 01013 //ignoreNonLocalEntries_ is true. 01014 if (_vec->Map().Comm().NumProc() < 2 || ignoreNonLocalEntries_) { 01015 return(0); 01016 } 01017 01018 01019 01020 //First build a map that describes the data in nonlocalIDs_/nonlocalCoefs_. 01021 //We'll use the arbitrary distribution constructor of Map. 01022 01023 Epetra_BlockMap sourceMap(-1, numNonlocalIDs_, 01024 nonlocalIDs_, nonlocalElementSize_, 01025 _vec->Map().IndexBase(), _vec->Map().Comm()); 01026 01027 //Now build a vector to hold our nonlocalCoefs_, and to act as the source- 01028 //vector for our import operation. 01029 Epetra_MultiVector nonlocalVector(sourceMap, 1); 01030 01031 int i,j; 01032 for(i=0; i<numNonlocalIDs_; ++i) { 01033 for(j=0; j<nonlocalElementSize_[i]; ++j) { 01034 nonlocalVector.ReplaceGlobalValue(nonlocalIDs_[i], j, 0, 01035 nonlocalCoefs_[i][j]); 01036 } 01037 } 01038 01039 Epetra_Export exporter(sourceMap, _vec->Map()); 01040 01041 EPETRA_CHK_ERR( _vec->Export(nonlocalVector, exporter, mode) ); 01042 01043 destroyNonlocalData(); 01044 01045 return(0); 01046 }
| void libMesh::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 163 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
00164 { 00165 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00166 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00167 00168 p.first++; 00169 }
| void libMesh::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 176 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
00177 { 00178 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00179 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00180 00181 p.second++; 00182 }
| virtual void libMesh::NumericVector< T >::init | ( | const NumericVector< T > & | other, | |
| const bool | fast = false | |||
| ) | [pure virtual, inherited] |
Creates a vector that has the same dimension and storage type as other, including ghost dofs.
| void libMesh::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.
Definition at line 732 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init(), libMesh::NumericVector< T >::local_size(), libMesh::NumericVector< T >::size(), and libMesh::NumericVector< T >::type().
00734 { 00735 this->init(other.size(),other.local_size(),fast,other.type()); 00736 }
| void libMesh::EpetraVector< T >::init | ( | const numeric_index_type | n, | |
| const numeric_index_type | n_local, | |||
| const std::vector< numeric_index_type > & | , | |||
| const bool | fast = false, |
|||
| const ParallelType | type = AUTOMATIC | |||
| ) | [inline, virtual] |
Create a vector that holds the local indices plus those specified in the ghost argument.
Implements libMesh::NumericVector< T >.
Definition at line 805 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
| void libMesh::EpetraVector< T >::init | ( | const numeric_index_type | N, | |
| const bool | fast = false, |
|||
| const ParallelType | type = AUTOMATIC | |||
| ) | [inline, virtual] |
call init with n_local = N,
Implements libMesh::NumericVector< T >.
Definition at line 819 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
| void libMesh::EpetraVector< T >::init | ( | const numeric_index_type | N, | |
| const numeric_index_type | 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 libMesh::NumericVector< T >.
Definition at line 751 of file trilinos_epetra_vector.h.
References libMesh::NumericVector< T >::_is_closed, libMesh::NumericVector< T >::_is_initialized, libMesh::EpetraVector< T >::_map, libMesh::NumericVector< T >::_type, libMesh::EpetraVector< T >::_vec, libMeshEnums::AUTOMATIC, libMesh::COMM_WORLD, libMeshEnums::GHOSTED, libMesh::EpetraVector< T >::myCoefs_, libMesh::EpetraVector< T >::myFirstID_, libMesh::EpetraVector< T >::myNumIDs_, libMeshEnums::PARALLEL, libMeshEnums::SERIAL, and libMesh::EpetraVector< T >::zero().
Referenced by libMesh::EpetraVector< T >::EpetraVector(), and libMesh::EpetraVector< T >::init().
00755 { 00756 // We default to allocating n_local local storage 00757 numeric_index_type my_n_local = n_local; 00758 00759 if (type == AUTOMATIC) 00760 { 00761 if (n == n_local) 00762 this->_type = SERIAL; 00763 else 00764 this->_type = PARALLEL; 00765 } 00766 else if (type == GHOSTED) 00767 { 00768 // We don't yet support GHOSTED Epetra vectors, so to get the 00769 // same functionality we need a SERIAL vector with local 00770 // storage allocated for every entry. 00771 this->_type = SERIAL; 00772 my_n_local = n; 00773 } 00774 else 00775 this->_type = type; 00776 00777 libmesh_assert ((this->_type==SERIAL && n==my_n_local) || 00778 this->_type==PARALLEL); 00779 00780 _map = new Epetra_Map(static_cast<int>(n), 00781 my_n_local, 00782 0, 00783 Epetra_MpiComm (libMesh::COMM_WORLD)); 00784 00785 _vec = new Epetra_Vector(*_map); 00786 00787 myFirstID_ = _vec->Map().MinMyGID(); 00788 myNumIDs_ = _vec->Map().NumMyElements(); 00789 00790 //Currently we impose the restriction that NumVectors==1, so we won't 00791 //need the LDA argument when calling ExtractView. Hence the "dummy" arg. 00792 int dummy; 00793 _vec->ExtractView(&myCoefs_, &dummy); 00794 00795 this->_is_initialized = true; 00796 this->_is_closed = true; 00797 00798 if (fast == false) 00799 this->zero (); 00800 }
| virtual bool libMesh::NumericVector< T >::initialized | ( | ) | const [inline, virtual, inherited] |
- Returns:
- true if the vector has been initialized, false otherwise.
Definition at line 111 of file numeric_vector.h.
Referenced by libMesh::PetscVector< T >::_get_array(), libMesh::PetscVector< T >::_restore_array(), libMesh::LaspackVector< T >::abs(), libMesh::DistributedVector< T >::abs(), libMesh::LaspackVector< T >::add(), libMesh::DistributedVector< T >::add(), libMesh::DistributedVector< T >::add_vector(), libMesh::ImplicitSystem::assemble(), libMesh::EpetraVector< T >::clear(), libMesh::PetscVector< T >::clear(), libMesh::LaspackVector< T >::clear(), libMesh::EpetraVector< T >::close(), libMesh::LaspackVector< T >::close(), libMesh::DistributedVector< T >::close(), libMesh::NumericVector< T >::compare(), libMesh::PetscVector< T >::create_subvector(), libMesh::LaspackVector< T >::dot(), libMesh::EpetraVector< T >::first_local_index(), libMesh::PetscVector< T >::first_local_index(), libMesh::LaspackVector< T >::first_local_index(), libMesh::DistributedVector< T >::first_local_index(), libMesh::NumericVector< T >::global_relative_compare(), libMesh::PetscVector< T >::init(), libMesh::DistributedVector< T >::init(), libMesh::DistributedVector< T >::insert(), libMesh::DistributedVector< T >::l1_norm(), libMesh::DistributedVector< T >::l2_norm(), libMesh::EpetraVector< T >::last_local_index(), libMesh::PetscVector< T >::last_local_index(), libMesh::LaspackVector< T >::last_local_index(), libMesh::DistributedVector< T >::last_local_index(), libMesh::DistributedVector< T >::linfty_norm(), libMesh::NumericVector< T >::local_relative_compare(), libMesh::EpetraVector< T >::local_size(), libMesh::PetscVector< T >::local_size(), libMesh::LaspackVector< T >::local_size(), libMesh::DistributedVector< T >::local_size(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::localize_to_one(), libMesh::PetscVector< T >::map_global_to_local_index(), libMesh::EpetraVector< T >::max(), libMesh::LaspackVector< T >::max(), libMesh::DistributedVector< T >::max(), libMesh::EpetraVector< T >::min(), libMesh::LaspackVector< T >::min(), libMesh::DistributedVector< T >::min(), libMesh::EpetraVector< T >::operator()(), libMesh::LaspackVector< T >::operator()(), libMesh::DistributedVector< T >::operator()(), libMesh::DistributedVector< T >::operator+=(), libMesh::DistributedVector< T >::operator-=(), libMesh::LaspackVector< T >::operator=(), libMesh::DistributedVector< T >::operator=(), libMesh::NumericVector< T >::print(), libMesh::NumericVector< T >::print_global(), libMesh::LaspackVector< T >::scale(), libMesh::DistributedVector< T >::scale(), libMesh::LaspackVector< T >::set(), libMesh::DistributedVector< T >::set(), libMesh::EpetraVector< T >::size(), libMesh::PetscVector< T >::size(), libMesh::LaspackVector< T >::size(), libMesh::DistributedVector< T >::size(), libMesh::DistributedVector< T >::sum(), libMesh::EpetraVector< T >::zero(), libMesh::LaspackVector< T >::zero(), and libMesh::DistributedVector< T >::zero().
00111 { return _is_initialized; }
| int libMesh::EpetraVector< T >::inputNonlocalValue | ( | int | GID, | |
| double | value, | |||
| bool | accumulate | |||
| ) | [inline, private] |
Definition at line 906 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::allocatedNonlocalLength_, libMesh::EpetraVector< T >::nonlocalCoefs_, libMesh::EpetraVector< T >::nonlocalElementSize_, libMesh::EpetraVector< T >::nonlocalIDs_, and libMesh::EpetraVector< T >::numNonlocalIDs_.
Referenced by libMesh::EpetraVector< T >::inputValues().
00907 { 00908 int insertPoint = -1; 00909 00910 //find offset of GID in nonlocalIDs_ 00911 int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_, 00912 insertPoint); 00913 if (offset >= 0) { 00914 //if offset >= 0 00915 // put value in nonlocalCoefs_[offset][0] 00916 00917 if (accumulate) { 00918 nonlocalCoefs_[offset][0] += value; 00919 } 00920 else { 00921 nonlocalCoefs_[offset][0] = value; 00922 } 00923 } 00924 else { 00925 //else 00926 // insert GID in nonlocalIDs_ 00927 // insert 1 in nonlocalElementSize_ 00928 // insert value in nonlocalCoefs_ 00929 00930 int tmp1 = numNonlocalIDs_; 00931 int tmp2 = allocatedNonlocalLength_; 00932 int tmp3 = allocatedNonlocalLength_; 00933 EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_, 00934 tmp1, tmp2) ); 00935 --tmp1; 00936 EPETRA_CHK_ERR( Epetra_Util_insert(1, insertPoint, nonlocalElementSize_, 00937 tmp1, tmp3) ); 00938 double* values = new double[1]; 00939 values[0] = value; 00940 EPETRA_CHK_ERR( Epetra_Util_insert(values, insertPoint, nonlocalCoefs_, 00941 numNonlocalIDs_, allocatedNonlocalLength_) ); 00942 } 00943 00944 return(0); 00945 }
| int libMesh::EpetraVector< T >::inputNonlocalValues | ( | int | GID, | |
| int | numValues, | |||
| const double * | values, | |||
| bool | accumulate | |||
| ) | [inline, private] |
Definition at line 949 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::allocatedNonlocalLength_, libMesh::EpetraVector< T >::nonlocalCoefs_, libMesh::EpetraVector< T >::nonlocalElementSize_, libMesh::EpetraVector< T >::nonlocalIDs_, and libMesh::EpetraVector< T >::numNonlocalIDs_.
Referenced by libMesh::EpetraVector< T >::inputValues().
00951 { 00952 int insertPoint = -1; 00953 00954 //find offset of GID in nonlocalIDs_ 00955 int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_, 00956 insertPoint); 00957 if (offset >= 0) { 00958 //if offset >= 0 00959 // put value in nonlocalCoefs_[offset][0] 00960 00961 if (numValues != nonlocalElementSize_[offset]) { 00962 cerr << "Epetra_FEVector ERROR: block-size for GID " << GID << " is " 00963 << numValues<<" which doesn't match previously set block-size of " 00964 << nonlocalElementSize_[offset] << endl; 00965 return(-1); 00966 } 00967 00968 if (accumulate) { 00969 for(int j=0; j<numValues; ++j) { 00970 nonlocalCoefs_[offset][j] += values[j]; 00971 } 00972 } 00973 else { 00974 for(int j=0; j<numValues; ++j) { 00975 nonlocalCoefs_[offset][j] = values[j]; 00976 } 00977 } 00978 } 00979 else { 00980 //else 00981 // insert GID in nonlocalIDs_ 00982 // insert numValues in nonlocalElementSize_ 00983 // insert values in nonlocalCoefs_ 00984 00985 int tmp1 = numNonlocalIDs_; 00986 int tmp2 = allocatedNonlocalLength_; 00987 int tmp3 = allocatedNonlocalLength_; 00988 EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_, 00989 tmp1, tmp2) ); 00990 --tmp1; 00991 EPETRA_CHK_ERR( Epetra_Util_insert(numValues, insertPoint, nonlocalElementSize_, 00992 tmp1, tmp3) ); 00993 double* newvalues = new double[numValues]; 00994 for(int j=0; j<numValues; ++j) { 00995 newvalues[j] = values[j]; 00996 } 00997 EPETRA_CHK_ERR( Epetra_Util_insert(newvalues, insertPoint, nonlocalCoefs_, 00998 numNonlocalIDs_, allocatedNonlocalLength_) ); 00999 } 01000 01001 return(0); 01002 }
| int libMesh::EpetraVector< T >::inputValues | ( | int | numIDs, | |
| const int * | GIDs, | |||
| const int * | numValuesPerID, | |||
| const double * | values, | |||
| bool | accumulate | |||
| ) | [inline, private] |
Definition at line 871 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::ignoreNonLocalEntries_, and libMesh::EpetraVector< T >::inputNonlocalValues().
00876 { 00877 int offset=0; 00878 for(int i=0; i<numIDs; ++i) { 00879 int numValues = numValuesPerID[i]; 00880 if (_vec->Map().MyGID(GIDs[i])) { 00881 if (accumulate) { 00882 for(int j=0; j<numValues; ++j) { 00883 _vec->SumIntoGlobalValue(GIDs[i], j, 0, values[offset+j]); 00884 } 00885 } 00886 else { 00887 for(int j=0; j<numValues; ++j) { 00888 _vec->ReplaceGlobalValue(GIDs[i], j, 0, values[offset+j]); 00889 } 00890 } 00891 } 00892 else { 00893 if (!ignoreNonLocalEntries_) { 00894 EPETRA_CHK_ERR( inputNonlocalValues(GIDs[i], numValues, 00895 &(values[offset]), accumulate) ); 00896 } 00897 } 00898 offset += numValues; 00899 } 00900 00901 return(0); 00902 }
| int libMesh::EpetraVector< T >::inputValues | ( | int | numIDs, | |
| const int * | GIDs, | |||
| const double * | values, | |||
| bool | accumulate | |||
| ) | [inline, private] |
Definition at line 842 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::ignoreNonLocalEntries_, and libMesh::EpetraVector< T >::inputNonlocalValue().
Referenced by libMesh::EpetraVector< T >::ReplaceGlobalValues(), and libMesh::EpetraVector< T >::SumIntoGlobalValues().
00846 { 00847 //Important note!! This method assumes that there is only 1 point 00848 //associated with each element. 00849 00850 for(int i=0; i<numIDs; ++i) { 00851 if (_vec->Map().MyGID(GIDs[i])) { 00852 if (accumulate) { 00853 _vec->SumIntoGlobalValue(GIDs[i], 0, 0, values[i]); 00854 } 00855 else { 00856 _vec->ReplaceGlobalValue(GIDs[i], 0, 0, values[i]); 00857 } 00858 } 00859 else { 00860 if (!ignoreNonLocalEntries_) { 00861 EPETRA_CHK_ERR( inputNonlocalValue(GIDs[i], values[i], accumulate) ); 00862 } 00863 } 00864 } 00865 00866 return(0); 00867 }
| virtual void libMesh::NumericVector< T >::insert | ( | const NumericVector< T > & | V, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [pure virtual, inherited] |
, where U and V are type NumericVector<T> and you want to specify WHERE to insert the NumericVector<T> V
| virtual void libMesh::EpetraVector< T >::insert | ( | const DenseSubVector< T > & | V, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [virtual] |
where V is type DenseSubVector<T> and you want to specify WHERE to insert it
Implements libMesh::NumericVector< T >.
| virtual void libMesh::EpetraVector< T >::insert | ( | const DenseVector< T > & | V, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [virtual] |
where V is type DenseVector<T> and you want to specify WHERE to insert it
Implements libMesh::NumericVector< T >.
| virtual void libMesh::EpetraVector< T >::insert | ( | const NumericVector< T > & | V, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [virtual] |
, where U and V are type NumericVector<T> and you want to specify WHERE to insert the NumericVector<T> V
| virtual void libMesh::EpetraVector< T >::insert | ( | const std::vector< T > & | v, | |
| const std::vector< numeric_index_type > & | dof_indices | |||
| ) | [virtual] |
where v is a DenseVector<T> and you want to specify WHERE to insert it
Implements libMesh::NumericVector< T >.
| Real libMesh::EpetraVector< T >::l1_norm | ( | ) | const [inline, virtual] |
- Returns:
- the
-norm of the vector, i.e. the sum of the absolute values.
Implements libMesh::NumericVector< T >.
Definition at line 67 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::NumericVector< T >::closed(), and libMesh::Real.
| Real libMesh::EpetraVector< T >::l2_norm | ( | ) | const [inline, virtual] |
- Returns:
- the
-norm of the vector, i.e. the square root of the sum of the squares of the elements.
Implements libMesh::NumericVector< T >.
Definition at line 79 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::NumericVector< T >::closed(), and libMesh::Real.
| numeric_index_type libMesh::EpetraVector< T >::last_local_index | ( | ) | const [inline, virtual] |
- Returns:
- the index of the last vector element actually stored on this processor
Implements libMesh::NumericVector< T >.
Definition at line 929 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, and libMesh::NumericVector< T >::initialized().
Referenced by libMesh::EpetraVector< T >::operator()().
00930 { 00931 libmesh_assert (this->initialized()); 00932 00933 return _vec->Map().MaxMyGID()+1; 00934 }
| Real libMesh::EpetraVector< T >::linfty_norm | ( | ) | const [inline, virtual] |
- Returns:
- the maximum absolute value of the elements of this vector, which is the
-norm of a vector.
Implements libMesh::NumericVector< T >.
Definition at line 91 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::NumericVector< T >::closed(), and libMesh::Real.
| int libMesh::NumericVector< T >::local_relative_compare | ( | const NumericVector< T > & | other_vector, | |
| const Real | threshold = TOLERANCE | |||
| ) | const [inline, virtual, inherited] |
- Returns:
-1whenthisis equivalent toother_vector, up to the given local relativethreshold. When differences occur, the return value contains the first index where the difference(a[i]-b[i])/max(a[i],b[i]) exceeded the threshold. When no threshold is given, thelibMeshTOLERANCEis used.
Definition at line 122 of file numeric_vector.C.
References std::abs(), libMesh::CommWorld, libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::last_local_index(), std::max(), libMesh::NumericVector< T >::max(), and libMesh::Parallel::Communicator::min().
00124 { 00125 libmesh_assert (this->initialized()); 00126 libmesh_assert (other_vector.initialized()); 00127 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00128 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00129 00130 int first_different_i = std::numeric_limits<int>::max(); 00131 numeric_index_type i = first_local_index(); 00132 00133 do 00134 { 00135 if ( std::abs( (*this)(i) - other_vector(i) ) > threshold * 00136 std::max(std::abs((*this)(i)), std::abs(other_vector(i)))) 00137 first_different_i = i; 00138 else 00139 i++; 00140 } 00141 while (first_different_i==std::numeric_limits<int>::max() 00142 && i<last_local_index()); 00143 00144 // Find the correct first differing index in parallel 00145 CommWorld.min(first_different_i); 00146 00147 if (first_different_i == std::numeric_limits<int>::max()) 00148 return -1; 00149 00150 return first_different_i; 00151 }
| numeric_index_type libMesh::EpetraVector< T >::local_size | ( | ) | const [inline, virtual] |
- Returns:
- the local size of the vector (index_stop-index_start)
Implements libMesh::NumericVector< T >.
Definition at line 909 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, and libMesh::NumericVector< T >::initialized().
Referenced by libMesh::EpetraVector< T >::localize(), libMesh::EpetraVector< T >::localize_to_one(), and libMesh::EpetraVector< T >::operator=().
00910 { 00911 libmesh_assert (this->initialized()); 00912 00913 return _vec->MyLength(); 00914 }
| virtual void libMesh::NumericVector< T >::localize | ( | NumericVector< T > & | v_local, | |
| const std::vector< numeric_index_type > & | send_list | |||
| ) | const [pure virtual, inherited] |
Creates a local vector v_local containing only information relevant to this processor, as defined by the send_list.
| virtual void libMesh::NumericVector< T >::localize | ( | NumericVector< T > & | v_local | ) | const [pure virtual, inherited] |
Same, but fills a NumericVector<T> instead of a std::vector.
| void libMesh::EpetraVector< T >::localize | ( | const numeric_index_type | first_local_idx, | |
| const numeric_index_type | last_local_idx, | |||
| const std::vector< numeric_index_type > & | send_list | |||
| ) | [virtual] |
Updates a local vector with selected values from neighboring processors, as defined by send_list.
Implements libMesh::NumericVector< T >.
| void libMesh::EpetraVector< T >::localize | ( | NumericVector< T > & | v_local, | |
| const std::vector< numeric_index_type > & | send_list | |||
| ) | const |
Creates a local vector v_local containing only information relevant to this processor, as defined by the send_list.
| void libMesh::EpetraVector< T >::localize | ( | NumericVector< T > & | v_local | ) | const [inline] |
Same, but fills a NumericVector<T> instead of a std::vector.
Definition at line 450 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_map, and libMesh::EpetraVector< T >::_vec.
00451 { 00452 EpetraVector<T>* v_local = libmesh_cast_ptr<EpetraVector<T>*>(&v_local_in); 00453 00454 Epetra_Map rootMap = Epetra_Util::Create_Root_Map( *_map, -1); 00455 v_local->_vec->ReplaceMap(rootMap); 00456 00457 Epetra_Import importer(v_local->_vec->Map(), *_map); 00458 v_local->_vec->Import(*_vec, importer, Insert); 00459 }
| void libMesh::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 libMesh::NumericVector< T >.
Definition at line 573 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::Parallel::Communicator::allgather(), libMesh::CommWorld, libMesh::EpetraVector< T >::local_size(), and libMesh::EpetraVector< T >::size().
00574 { 00575 // This function must be run on all processors at once 00576 parallel_only(); 00577 00578 const unsigned int n = this->size(); 00579 const unsigned int nl = this->local_size(); 00580 00581 libmesh_assert(this->_vec); 00582 00583 v_local.clear(); 00584 v_local.reserve(n); 00585 00586 // build up my local part 00587 for (unsigned int i=0; i<nl; i++) 00588 v_local.push_back((*this->_vec)[i]); 00589 00590 CommWorld.allgather (v_local); 00591 }
| void libMesh::EpetraVector< T >::localize_to_one | ( | std::vector< T > & | v_local, | |
| const processor_id_type | 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 libMesh::NumericVector< T >.
Definition at line 596 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::CommWorld, libMesh::Parallel::Communicator::gather(), libMesh::EpetraVector< T >::local_size(), libMesh::n_processors(), and libMesh::EpetraVector< T >::size().
00598 { 00599 // This function must be run on all processors at once 00600 parallel_only(); 00601 00602 const unsigned int n = this->size(); 00603 const unsigned int nl = this->local_size(); 00604 00605 libmesh_assert_less (pid, libMesh::n_processors()); 00606 libmesh_assert(this->_vec); 00607 00608 v_local.clear(); 00609 v_local.reserve(n); 00610 00611 00612 // build up my local part 00613 for (unsigned int i=0; i<nl; i++) 00614 v_local.push_back((*this->_vec)[i]); 00615 00616 CommWorld.gather (pid, v_local); 00617 }
| Real libMesh::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 libMesh::NumericVector< T >.
Definition at line 967 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, and libMesh::NumericVector< T >::initialized().
00968 { 00969 libmesh_assert (this->initialized()); 00970 00971 T value; 00972 00973 _vec->MaxValue(&value); 00974 00975 return value; 00976 }
| Real libMesh::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 libMesh::NumericVector< T >.
Definition at line 952 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, and libMesh::NumericVector< T >::initialized().
00953 { 00954 libmesh_assert (this->initialized()); 00955 00956 T value; 00957 00958 _vec->MinValue(&value); 00959 00960 return value; 00961 }
| static unsigned int libMesh::ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 79 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
00080 { return _n_objects; }
| T libMesh::EpetraVector< T >::operator() | ( | const numeric_index_type | i | ) | const [inline, virtual] |
Access components, returns U(i).
Implements libMesh::NumericVector< T >.
Definition at line 939 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), and libMesh::EpetraVector< T >::last_local_index().
00940 { 00941 libmesh_assert (this->initialized()); 00942 libmesh_assert ( ((i >= this->first_local_index()) && 00943 (i < this->last_local_index())) ); 00944 00945 return (*_vec)[i-this->first_local_index()]; 00946 }
| NumericVector<T>& libMesh::NumericVector< T >::operator*= | ( | const T | a | ) | [inline, inherited] |
Multiplication operator. Equivalent to U.scale(a)
Definition at line 350 of file numeric_vector.h.
00350 { this->scale(a); return *this; }
| virtual NumericVector<T>& libMesh::NumericVector< T >::operator+= | ( | const NumericVector< T > & | V | ) | [pure virtual, inherited] |
Addition operator. Fast equivalent to U.add(1, V).
| NumericVector< T > & libMesh::EpetraVector< T >::operator+= | ( | const NumericVector< T > & | V | ) | [inline] |
Addition operator. Fast equivalent to U.add(1, V).
Definition at line 104 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::add(), and libMesh::NumericVector< T >::closed().
| virtual NumericVector<T>& libMesh::NumericVector< T >::operator-= | ( | const NumericVector< T > & | V | ) | [pure virtual, inherited] |
Subtraction operator. Fast equivalent to U.add(-1, V).
| NumericVector< T > & libMesh::EpetraVector< T >::operator-= | ( | const NumericVector< T > & | V | ) | [inline] |
Subtraction operator. Fast equivalent to U.add(-1, V).
Definition at line 117 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::add(), and libMesh::NumericVector< T >::closed().
| NumericVector<T>& libMesh::NumericVector< T >::operator/= | ( | const T | a | ) | [inline, inherited] |
Division operator. Equivalent to U.scale(1./a)
Definition at line 356 of file numeric_vector.h.
00356 { this->scale(1./a); return *this; }
| NumericVector< T > & libMesh::EpetraVector< T >::operator= | ( | const std::vector< T > & | v | ) | [inline, virtual] |
: 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 libMesh::NumericVector< T >.
Definition at line 413 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::first_local_index(), libMesh::EpetraVector< T >::local_size(), and libMesh::EpetraVector< T >::size().
00414 { 00415 T * values = _vec->Values(); 00416 00421 if(this->size() == v.size()) 00422 { 00423 const unsigned int nl=this->local_size(); 00424 const unsigned int fli=this->first_local_index(); 00425 00426 for(unsigned int i=0;i<nl;i++) 00427 values[i]=v[fli+i]; 00428 } 00429 00434 else 00435 { 00436 libmesh_assert_equal_to (v.size(), this->local_size()); 00437 00438 const unsigned int nl=this->local_size(); 00439 00440 for(unsigned int i=0;i<nl;i++) 00441 values[i]=v[i]; 00442 } 00443 00444 return *this; 00445 }
| EpetraVector< T > & libMesh::EpetraVector< T >::operator= | ( | const EpetraVector< T > & | V | ) | [inline, virtual] |
: copy all components.
Implements libMesh::NumericVector< T >.
Definition at line 402 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
| NumericVector< T > & libMesh::EpetraVector< T >::operator= | ( | const NumericVector< T > & | V | ) | [inline] |
: copy all components.
Definition at line 389 of file trilinos_epetra_vector.C.
| NumericVector< T > & libMesh::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).
: fill all components.
Implements libMesh::NumericVector< T >.
Definition at line 378 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
00379 { 00380 _vec->PutScalar(s_in); 00381 00382 return *this; 00383 }
| virtual void libMesh::NumericVector< T >::pointwise_mult | ( | const NumericVector< T > & | vec1, | |
| const NumericVector< T > & | vec2 | |||
| ) | [pure virtual, inherited] |
Computes the pointwise (i.e. component-wise) product of vec1 and vec2 and stores the result in *this.
Referenced by libMesh::TensorShellMatrix< T >::get_diagonal().
| void libMesh::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.
Definition at line 366 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
00368 { 00369 const EpetraVector<T>* V1 = libmesh_cast_ptr<const EpetraVector<T>*>(&vec1); 00370 const EpetraVector<T>* V2 = libmesh_cast_ptr<const EpetraVector<T>*>(&vec2); 00371 00372 _vec->Multiply(1.0, *V1->_vec, *V2->_vec, 0.0); 00373 }
| void libMesh::NumericVector< Complex >::print | ( | std::ostream & | os | ) | const [inline, inherited] |
Definition at line 777 of file numeric_vector.h.
References libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::NumericVector< T >::local_size(), and libMesh::NumericVector< T >::size().
00778 { 00779 libmesh_assert (this->initialized()); 00780 os << "Size\tglobal = " << this->size() 00781 << "\t\tlocal = " << this->local_size() << std::endl; 00782 00783 // std::complex<>::operator<<() is defined, but use this form 00784 os << "#\tReal part\t\tImaginary part" << std::endl; 00785 for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++) 00786 os << i << "\t" 00787 << (*this)(i).real() << "\t\t" 00788 << (*this)(i).imag() << std::endl; 00789 }
| void libMesh::NumericVector< T >::print | ( | std::ostream & | os = libMesh::out |
) | const [inline, virtual, inherited] |
Prints the local contents of the vector, by default to libMesh::out
Definition at line 795 of file numeric_vector.h.
References libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::NumericVector< T >::local_size(), and libMesh::NumericVector< T >::size().
00796 { 00797 libmesh_assert (this->initialized()); 00798 os << "Size\tglobal = " << this->size() 00799 << "\t\tlocal = " << this->local_size() << std::endl; 00800 00801 os << "#\tValue" << std::endl; 00802 for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++) 00803 os << i << "\t" << (*this)(i) << std::endl; 00804 }
| void libMesh::NumericVector< Complex >::print_global | ( | std::ostream & | os | ) | const [inline, inherited] |
Definition at line 810 of file numeric_vector.h.
References libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::localize(), libMesh::processor_id(), and libMesh::NumericVector< T >::size().
00811 { 00812 libmesh_assert (this->initialized()); 00813 00814 std::vector<Complex> v(this->size()); 00815 this->localize(v); 00816 00817 // Right now we only want one copy of the output 00818 if (libMesh::processor_id()) 00819 return; 00820 00821 os << "Size\tglobal = " << this->size() << std::endl; 00822 os << "#\tReal part\t\tImaginary part" << std::endl; 00823 for (numeric_index_type i=0; i!=v.size(); i++) 00824 os << i << "\t" 00825 << v[i].real() << "\t\t" 00826 << v[i].imag() << std::endl; 00827 }
| void libMesh::NumericVector< T >::print_global | ( | std::ostream & | os = libMesh::out |
) | const [inline, virtual, inherited] |
Prints the global contents of the vector, by default to libMesh::out
Definition at line 832 of file numeric_vector.h.
References libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::localize(), libMesh::processor_id(), and libMesh::NumericVector< T >::size().
00833 { 00834 libmesh_assert (this->initialized()); 00835 00836 std::vector<T> v(this->size()); 00837 this->localize(v); 00838 00839 // Right now we only want one copy of the output 00840 if (libMesh::processor_id()) 00841 return; 00842 00843 os << "Size\tglobal = " << this->size() << std::endl; 00844 os << "#\tValue" << std::endl; 00845 for (numeric_index_type i=0; i!=v.size(); i++) 00846 os << i << "\t" << v[i] << std::endl; 00847 }
| void libMesh::ReferenceCounter::print_info | ( | std::ostream & | out = libMesh::out |
) | [static, inherited] |
Prints the reference information, by default to libMesh::out.
Definition at line 88 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
00089 { 00090 if( _enable_print_counter ) out_stream << ReferenceCounter::get_info(); 00091 }
| void libMesh::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 libMesh::NumericVector< T >.
Definition at line 622 of file trilinos_epetra_vector.C.
00623 { 00624 libmesh_not_implemented(); 00625 00626 // libmesh_assert (this->initialized()); 00627 // libmesh_assert (this->closed()); 00628 00629 // int ierr=0; 00630 // EpetraViewer epetra_viewer; 00631 00632 00633 // ierr = EpetraViewerCreate (libMesh::COMM_WORLD, 00634 // &epetra_viewer); 00635 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00636 00637 // /** 00638 // * Create an ASCII file containing the matrix 00639 // * if a filename was provided. 00640 // */ 00641 // if (name != "NULL") 00642 // { 00643 // ierr = EpetraViewerASCIIOpen( libMesh::COMM_WORLD, 00644 // name.c_str(), 00645 // &epetra_viewer); 00646 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00647 00648 // ierr = EpetraViewerSetFormat (epetra_viewer, 00649 // EPETRA_VIEWER_ASCII_MATLAB); 00650 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00651 00652 // ierr = VecView (_vec, epetra_viewer); 00653 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00654 // } 00655 00656 // /** 00657 // * Otherwise the matrix will be dumped to the screen. 00658 // */ 00659 // else 00660 // { 00661 // ierr = EpetraViewerSetFormat (EPETRA_VIEWER_STDOUT_WORLD, 00662 // EPETRA_VIEWER_ASCII_MATLAB); 00663 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00664 00665 // ierr = VecView (_vec, EPETRA_VIEWER_STDOUT_WORLD); 00666 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00667 // } 00668 00669 00670 // /** 00671 // * Destroy the viewer. 00672 // */ 00673 // ierr = EpetraViewerDestroy (epetra_viewer); 00674 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00675 }
| void libMesh::EpetraVector< T >::reciprocal | ( | ) | [inline, virtual] |
Replace each entry v_i of this vector by its reciprocal, 1/v_i.
Implements libMesh::NumericVector< T >.
Definition at line 144 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, std::abs(), libMesh::EpetraVector< T >::close(), libMesh::err, and std::min().
00145 { 00146 // The Epetra::reciprocal() function takes a constant reference to *another* vector, 00147 // and fills _vec with its reciprocal. Does that mean we can't pass *_vec as the 00148 // argument? 00149 // _vec->reciprocal( *_vec ); 00150 00151 // Alternatively, compute the reciprocal by hand... see also the add(T) member that does this... 00152 const unsigned int nl = _vec->MyLength(); 00153 00154 T* values = _vec->Values(); 00155 00156 for (unsigned int i=0; i<nl; i++) 00157 { 00158 // Don't divide by zero (maybe only check this in debug mode?) 00159 if (std::abs(values[i]) < std::numeric_limits<T>::min()) 00160 { 00161 libMesh::err << "Error, divide by zero in DistributedVector<T>::reciprocal()!" << std::endl; 00162 libmesh_error(); 00163 } 00164 00165 values[i] = 1. / values[i]; 00166 } 00167 00168 // Leave the vector in a closed state... 00169 this->close(); 00170 }
| int libMesh::EpetraVector< T >::ReplaceGlobalValues | ( | int | numIDs, | |
| const int * | GIDs, | |||
| const int * | numValuesPerID, | |||
| const double * | values | |||
| ) | [inline, private] |
Definition at line 833 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::inputValues().
00836 { 00837 return( inputValues( numIDs, GIDs, numValuesPerID, values, false) ); 00838 }
| int libMesh::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 821 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::inputValues().
00823 { 00824 if (GIDs.Length() != values.Length()) { 00825 return(-1); 00826 } 00827 00828 return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), false) ); 00829 }
| int libMesh::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 813 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::inputValues().
00815 { 00816 return( inputValues( numIDs, GIDs, values, false) ); 00817 }
| void libMesh::EpetraVector< T >::scale | ( | const T | factor | ) | [inline, virtual] |
Scale each element of the vector by the given factor.
Implements libMesh::NumericVector< T >.
Definition at line 340 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
00341 { 00342 _vec->Scale(factor_in); 00343 }
| void libMesh::EpetraVector< T >::set | ( | const numeric_index_type | i, | |
| const T | value | |||
| ) | [virtual] |
v(i) = value
Implements libMesh::NumericVector< T >.
| void libMesh::EpetraVector< T >::setIgnoreNonLocalEntries | ( | bool | flag | ) | [inline, private] |
Set whether or not non-local data values should be ignored.
Definition at line 573 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::ignoreNonLocalEntries_.
00573 { 00574 ignoreNonLocalEntries_ = flag; 00575 }
| numeric_index_type libMesh::EpetraVector< T >::size | ( | ) | const [inline, virtual] |
- Returns:
- dimension of the vector. This function was formerly called
n(), but was renamed to get theEpetraVector<T>class closer to the C++ standard library'sstd::vectorcontainer.
Implements libMesh::NumericVector< T >.
Definition at line 898 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, and libMesh::NumericVector< T >::initialized().
Referenced by libMesh::EpetraVector< T >::add(), libMesh::EpetraVector< T >::localize(), libMesh::EpetraVector< T >::localize_to_one(), and libMesh::EpetraVector< T >::operator=().
00899 { 00900 libmesh_assert (this->initialized()); 00901 00902 return _vec->GlobalLength(); 00903 }
| Real libMesh::NumericVector< T >::subset_l1_norm | ( | const std::set< numeric_index_type > & | indices | ) | const [inline, virtual, inherited] |
- Returns:
- the
-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 298 of file numeric_vector.C.
References std::abs(), libMesh::CommWorld, libMesh::Real, and libMesh::Parallel::Communicator::sum().
Referenced by libMesh::System::discrete_var_norm().
00299 { 00300 const NumericVector<T> & v = *this; 00301 00302 std::set<numeric_index_type>::const_iterator it = indices.begin(); 00303 const std::set<numeric_index_type>::const_iterator it_end = indices.end(); 00304 00305 Real norm = 0; 00306 00307 for(; it!=it_end; ++it) 00308 norm += std::abs(v(*it)); 00309 00310 CommWorld.sum(norm); 00311 00312 return norm; 00313 }
| Real libMesh::NumericVector< T >::subset_l2_norm | ( | const std::set< numeric_index_type > & | indices | ) | const [inline, virtual, inherited] |
- Returns:
- the
-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 316 of file numeric_vector.C.
References libMesh::CommWorld, libMesh::TensorTools::norm_sq(), libMesh::Real, and libMesh::Parallel::Communicator::sum().
Referenced by libMesh::System::discrete_var_norm().
00317 { 00318 const NumericVector<T> & v = *this; 00319 00320 std::set<numeric_index_type>::const_iterator it = indices.begin(); 00321 const std::set<numeric_index_type>::const_iterator it_end = indices.end(); 00322 00323 Real norm = 0; 00324 00325 for(; it!=it_end; ++it) 00326 norm += TensorTools::norm_sq(v(*it)); 00327 00328 CommWorld.sum(norm); 00329 00330 return std::sqrt(norm); 00331 }
| Real libMesh::NumericVector< T >::subset_linfty_norm | ( | const std::set< numeric_index_type > & | indices | ) | const [inline, virtual, inherited] |
- Returns:
- the maximum absolute value of the specified entries of this vector, which is the
-norm of a vector.
Note that the indices must necessary live on this processor.
Definition at line 334 of file numeric_vector.C.
References libMesh::NumericVector< T >::abs(), libMesh::CommWorld, libMesh::Parallel::Communicator::max(), and libMesh::Real.
Referenced by libMesh::System::discrete_var_norm().
00335 { 00336 const NumericVector<T> & v = *this; 00337 00338 std::set<numeric_index_type>::const_iterator it = indices.begin(); 00339 const std::set<numeric_index_type>::const_iterator it_end = indices.end(); 00340 00341 Real norm = 0; 00342 00343 for(; it!=it_end; ++it) 00344 { 00345 Real value = std::abs(v(*it)); 00346 if(value > norm) 00347 norm = value; 00348 } 00349 00350 CommWorld.max(norm); 00351 00352 return norm; 00353 }
| T libMesh::EpetraVector< T >::sum | ( | ) | const [inline, virtual] |
- Returns:
- the sum of values in a vector
Implements libMesh::NumericVector< T >.
Definition at line 48 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::NumericVector< T >::closed(), libMesh::CommWorld, and libMesh::Parallel::Communicator::sum().
00049 { 00050 libmesh_assert(this->closed()); 00051 00052 const unsigned int nl = _vec->MyLength(); 00053 00054 T sum=0.0; 00055 00056 T * values = _vec->Values(); 00057 00058 for (unsigned int i=0; i<nl; i++) 00059 sum += values[i]; 00060 00061 CommWorld.sum<T>(sum); 00062 00063 return sum; 00064 }
| int libMesh::EpetraVector< T >::SumIntoGlobalValues | ( | int | numIDs, | |
| const int * | GIDs, | |||
| const int * | numValuesPerID, | |||
| const double * | values | |||
| ) | [inline, private] |
Definition at line 804 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::inputValues().
00807 { 00808 return( inputValues( numIDs, GIDs, numValuesPerID, values, true) ); 00809 }
| int libMesh::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 792 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::inputValues().
00794 { 00795 if (GIDs.Length() != values.Length()) { 00796 return(-1); 00797 } 00798 00799 return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), true) ); 00800 }
| int libMesh::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 784 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::inputValues().
00786 { 00787 return( inputValues( numIDs, GIDs, values, true) ); 00788 }
| void libMesh::NumericVector< T >::swap | ( | NumericVector< T > & | v | ) | [inline, virtual, inherited] |
Exchanges the values/sizes of two vectors. There should be enough indirection in subclasses to make this an O(1) header-swap operation.
Definition at line 853 of file numeric_vector.h.
References libMesh::NumericVector< T >::_is_closed, libMesh::NumericVector< T >::_is_initialized, and libMesh::NumericVector< T >::_type.
Referenced by libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), and libMesh::AdjointRefinementEstimator::estimate_error().
00854 { 00855 std::swap(_is_closed, v._is_closed); 00856 std::swap(_is_initialized, v._is_initialized); 00857 std::swap(_type, v._type); 00858 }
| void libMesh::EpetraVector< T >::swap | ( | NumericVector< T > & | v | ) | [inline, virtual] |
Swaps the raw Epetra vector context pointers.
Definition at line 982 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_destroy_vec_on_exit, and libMesh::EpetraVector< T >::_vec.
Referenced by libMesh::Problem_Interface::computeF().
00983 { 00984 NumericVector<T>::swap(other); 00985 00986 EpetraVector<T>& v = libmesh_cast_ref<EpetraVector<T>&>(other); 00987 00988 std::swap(_vec, v._vec); 00989 std::swap(_destroy_vec_on_exit, v._destroy_vec_on_exit); 00990 }
| ParallelType& libMesh::NumericVector< T >::type | ( | ) | [inline, inherited] |
- Returns:
- the type (SERIAL, PARALLEL, GHOSTED) of the vector.
Definition at line 121 of file numeric_vector.h.
00121 { return _type; }
| ParallelType libMesh::NumericVector< T >::type | ( | ) | const [inline, inherited] |
- Returns:
- the type (SERIAL, PARALLEL, GHOSTED) of the vector.
Definition at line 116 of file numeric_vector.h.
Referenced by libMesh::PetscVector< T >::_get_array(), libMesh::PetscVector< T >::_restore_array(), libMesh::PetscVector< T >::abs(), libMesh::PetscVector< T >::add(), libMesh::PetscVector< T >::close(), libMesh::DofMap::enforce_constraints_exactly(), libMesh::PetscVector< T >::get(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), libMesh::EpetraVector< T >::init(), libMesh::LaspackVector< T >::init(), libMesh::DistributedVector< T >::init(), libMesh::PetscVector< T >::localize(), libMesh::PetscVector< T >::operator()(), libMesh::MeshFunction::operator()(), libMesh::PetscVector< T >::operator=(), libMesh::PetscVector< T >::pointwise_mult(), libMesh::System::project_vector(), libMesh::System::read_serialized_vector(), libMesh::PetscVector< T >::scale(), and libMesh::PetscVector< T >::zero().
00116 { return _type; }
| Epetra_Vector* libMesh::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 LibMeshVecDestroy()!
Definition at line 490 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec.
Referenced by libMesh::EpetraMatrix< T >::get_diagonal(), and libMesh::NoxNonlinearSolver< T >::solve().
| void libMesh::EpetraVector< T >::zero | ( | ) | [inline, virtual] |
Set all entries to zero. Equivalent to v = 0, but more obvious and faster.
Implements libMesh::NumericVector< T >.
Definition at line 858 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, libMesh::NumericVector< T >::closed(), and libMesh::NumericVector< T >::initialized().
Referenced by libMesh::Problem_Interface::computeF(), and libMesh::EpetraVector< T >::init().
00859 { 00860 libmesh_assert (this->initialized()); 00861 libmesh_assert (this->closed()); 00862 00863 _vec->PutScalar(0.0); 00864 }
| AutoPtr< NumericVector< T > > libMesh::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 libMesh::NumericVector< T >.
Definition at line 870 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::add_vector().
Friends And Related Function Documentation
| 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 583 of file numeric_vector.h.
Member Data Documentation
ReferenceCounter::Counts libMesh::ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 118 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
bool libMesh::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 509 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::clear(), and libMesh::EpetraVector< T >::swap().
bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited] |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 137 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
bool libMesh::NumericVector< T >::_is_closed [protected, inherited] |
Flag to see if the Numeric assemble routines have been called yet
Definition at line 628 of file numeric_vector.h.
Referenced by libMesh::EpetraVector< T >::add(), libMesh::PetscVector< T >::add(), libMesh::LaspackVector< T >::add(), libMesh::EpetraVector< T >::clear(), libMesh::PetscVector< T >::clear(), libMesh::NumericVector< T >::clear(), libMesh::LaspackVector< T >::clear(), libMesh::DistributedVector< T >::clear(), libMesh::EpetraVector< T >::close(), libMesh::PetscVector< T >::close(), libMesh::LaspackVector< T >::close(), libMesh::DistributedVector< T >::close(), libMesh::NumericVector< Number >::closed(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< T >::init(), libMesh::DistributedVector< T >::localize(), libMesh::LaspackVector< T >::operator=(), libMesh::DistributedVector< T >::operator=(), libMesh::PetscVector< T >::PetscVector(), libMesh::PetscVector< T >::set(), libMesh::LaspackVector< T >::set(), and libMesh::NumericVector< T >::swap().
bool libMesh::NumericVector< T >::_is_initialized [protected, inherited] |
Flag to tell if init has been called yet
Definition at line 634 of file numeric_vector.h.
Referenced by libMesh::EpetraVector< T >::clear(), libMesh::PetscVector< T >::clear(), libMesh::NumericVector< T >::clear(), libMesh::LaspackVector< T >::clear(), libMesh::DistributedVector< T >::clear(), libMesh::PetscVector< T >::create_subvector(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< T >::init(), libMesh::DistributedVector< T >::init(), libMesh::NumericVector< Number >::initialized(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::operator=(), libMesh::PetscVector< T >::PetscVector(), and libMesh::NumericVector< T >::swap().
Epetra_Map* libMesh::EpetraVector< T >::_map [private] |
Holds the distributed Map
Definition at line 503 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::init(), and libMesh::EpetraVector< T >::localize().
Threads::spin_mutex libMesh::ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 131 of file reference_counter.h.
Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 126 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
ParallelType libMesh::NumericVector< T >::_type [protected, inherited] |
Type of vector
Definition at line 639 of file numeric_vector.h.
Referenced by libMesh::DistributedVector< T >::DistributedVector(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< T >::init(), libMesh::DistributedVector< T >::init(), libMesh::LaspackVector< T >::LaspackVector(), libMesh::PetscVector< T >::operator=(), libMesh::PetscVector< T >::PetscVector(), libMesh::NumericVector< T >::swap(), and libMesh::NumericVector< Number >::type().
Epetra_Vector* libMesh::EpetraVector< T >::_vec [private] |
Actual Epetra vector datatype to hold vector entries
Definition at line 498 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::abs(), libMesh::EpetraVector< T >::add(), libMesh::EpetraVector< T >::add_vector(), libMesh::EpetraVector< T >::clear(), libMesh::EpetraVector< T >::dot(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::EpetraVector< T >::FEoperatorequals(), libMesh::EpetraVector< T >::first_local_index(), libMesh::EpetraVector< T >::GlobalAssemble(), libMesh::EpetraVector< T >::init(), libMesh::EpetraVector< T >::inputValues(), libMesh::EpetraVector< T >::l1_norm(), libMesh::EpetraVector< T >::l2_norm(), libMesh::EpetraVector< T >::last_local_index(), libMesh::EpetraVector< T >::linfty_norm(), libMesh::EpetraVector< T >::local_size(), libMesh::EpetraVector< T >::localize(), libMesh::EpetraVector< T >::localize_to_one(), libMesh::EpetraVector< T >::max(), libMesh::EpetraVector< T >::min(), libMesh::EpetraVector< T >::operator()(), libMesh::EpetraVector< T >::operator=(), libMesh::EpetraVector< T >::pointwise_mult(), libMesh::EpetraVector< T >::reciprocal(), libMesh::EpetraVector< T >::scale(), libMesh::EpetraVector< T >::size(), libMesh::EpetraVector< T >::sum(), libMesh::EpetraVector< T >::swap(), libMesh::EpetraVector< T >::vec(), and libMesh::EpetraVector< T >::zero().
int libMesh::EpetraVector< T >::allocatedNonlocalLength_ [private] |
Definition at line 602 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::destroyNonlocalData(), libMesh::EpetraVector< T >::FEoperatorequals(), libMesh::EpetraVector< T >::inputNonlocalValue(), and libMesh::EpetraVector< T >::inputNonlocalValues().
bool libMesh::EpetraVector< T >::ignoreNonLocalEntries_ [private] |
Definition at line 605 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::GlobalAssemble(), libMesh::EpetraVector< T >::inputValues(), and libMesh::EpetraVector< T >::setIgnoreNonLocalEntries().
double* libMesh::EpetraVector< T >::myCoefs_ [private] |
Definition at line 597 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::EpetraVector(), and libMesh::EpetraVector< T >::init().
int libMesh::EpetraVector< T >::myFirstID_ [private] |
Definition at line 595 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::EpetraVector(), and libMesh::EpetraVector< T >::init().
int libMesh::EpetraVector< T >::myNumIDs_ [private] |
Definition at line 596 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::EpetraVector(), and libMesh::EpetraVector< T >::init().
double** libMesh::EpetraVector< T >::nonlocalCoefs_ [private] |
Definition at line 603 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::destroyNonlocalData(), libMesh::EpetraVector< T >::FEoperatorequals(), libMesh::EpetraVector< T >::GlobalAssemble(), libMesh::EpetraVector< T >::inputNonlocalValue(), and libMesh::EpetraVector< T >::inputNonlocalValues().
int* libMesh::EpetraVector< T >::nonlocalElementSize_ [private] |
Definition at line 600 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::destroyNonlocalData(), libMesh::EpetraVector< T >::FEoperatorequals(), libMesh::EpetraVector< T >::GlobalAssemble(), libMesh::EpetraVector< T >::inputNonlocalValue(), and libMesh::EpetraVector< T >::inputNonlocalValues().
int* libMesh::EpetraVector< T >::nonlocalIDs_ [private] |
Definition at line 599 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::destroyNonlocalData(), libMesh::EpetraVector< T >::FEoperatorequals(), libMesh::EpetraVector< T >::GlobalAssemble(), libMesh::EpetraVector< T >::inputNonlocalValue(), and libMesh::EpetraVector< T >::inputNonlocalValues().
int libMesh::EpetraVector< T >::numNonlocalIDs_ [private] |
Definition at line 601 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::destroyNonlocalData(), libMesh::EpetraVector< T >::FEoperatorequals(), libMesh::EpetraVector< T >::GlobalAssemble(), libMesh::EpetraVector< T >::inputNonlocalValue(), and libMesh::EpetraVector< T >::inputNonlocalValues().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:42 UTC
Hosted By: