$extrastylesheet
#include <trilinos_epetra_vector.h>

Public Member Functions | |
| EpetraVector (const Parallel::Communicator &comm, const ParallelType type=AUTOMATIC) | |
| EpetraVector (const Parallel::Communicator &comm, const numeric_index_type n, const ParallelType type=AUTOMATIC) | |
| EpetraVector (const Parallel::Communicator &comm, const numeric_index_type n, const numeric_index_type n_local, const ParallelType type=AUTOMATIC) | |
| EpetraVector (const Parallel::Communicator &comm, 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, const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD) | |
| ~EpetraVector () | |
| void | close () |
| void | clear () |
| void | zero () |
| virtual UniquePtr < NumericVector< T > > | zero_clone () const |
| UniquePtr< 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 NumericVector< T > & | operator/= (NumericVector< T > &v) |
| virtual void | reciprocal () |
| virtual void | conjugate () |
| 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 T *v, const std::vector< numeric_index_type > &dof_indices) |
| void | add_vector (const NumericVector< T > &V, const SparseMatrix< T > &A) |
| void | add_vector_transpose (const NumericVector< T > &V, const SparseMatrix< T > &A_trans) |
| virtual void | insert (const 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) |
| 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 Real | subset_l1_norm (const std::set< numeric_index_type > &indices) const |
| virtual Real | subset_l2_norm (const std::set< numeric_index_type > &indices) const |
| virtual Real | subset_linfty_norm (const std::set< numeric_index_type > &indices) const |
| virtual T | el (const numeric_index_type i) const |
| virtual void | get (const std::vector< numeric_index_type > &index, T *values) const |
| void | get (const std::vector< numeric_index_type > &index, std::vector< T > &values) const |
| NumericVector< T > & | operator*= (const T a) |
| NumericVector< T > & | operator/= (const T a) |
| void | add_vector (const std::vector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
| virtual void | add_vector (const NumericVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| void | add_vector (const DenseVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| void | add_vector (const NumericVector< T > &v, const ShellMatrix< T > &a) |
| void | insert (const std::vector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
| virtual void | insert (const NumericVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| void | insert (const DenseVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| void | insert (const DenseSubVector< T > &V, const std::vector< numeric_index_type > &dof_indices) |
| virtual int | compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
| virtual int | local_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
| virtual int | global_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
| virtual void | print (std::ostream &os=libMesh::out) const |
| template<> | |
| void | print (std::ostream &os) const |
| virtual void | print_global (std::ostream &os=libMesh::out) const |
| template<> | |
| void | print_global (std::ostream &os) const |
| virtual void | print_matlab (const std::string &="") const |
| const Parallel::Communicator & | comm () const |
| processor_id_type | n_processors () const |
| processor_id_type | processor_id () const |
Static Public Member Functions | |
| static UniquePtr < NumericVector< T > > | build (const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package()) |
| static UniquePtr < NumericVector< T > > | build (const SolverPackage solver_package=libMesh::default_solver_package()) |
| static std::string | get_info () |
| static void | print_info (std::ostream &out=libMesh::out) |
| static unsigned int | n_objects () |
| static void | enable_print_counter_info () |
| static void | disable_print_counter_info () |
Protected Types | |
| typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
| void | increment_constructor_count (const std::string &name) |
| void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
| bool | _is_closed |
| bool | _is_initialized |
| ParallelType | _type |
| const Parallel::Communicator & | _communicator |
Static Protected Attributes | |
| static Counts | _counts |
| static Threads::atomic < unsigned int > | _n_objects |
| static Threads::spin_mutex | _mutex |
| static bool | _enable_print_counter = true |
Private 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_ |
| unsigned char | last_edit |
| bool | ignoreNonLocalEntries_ |
Friends | |
| std::ostream & | operator<< (std::ostream &os, const NumericVector< T > &v) |
Epetra vector. Provides a nice interface to the Trilinos Epetra data structures for parallel vectors.
Definition at line 60 of file trilinos_epetra_vector.h.
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.
| libMesh::EpetraVector< T >::EpetraVector | ( | const Parallel::Communicator & | comm, |
| const ParallelType | type = AUTOMATIC |
||
| ) | [inline, explicit] |
Dummy-Constructor. Dimension=0
Definition at line 600 of file trilinos_epetra_vector.h.
References libMesh::NumericVector< T >::_type, and libMesh::NumericVector< T >::type().
: NumericVector<T>(comm, type), _destroy_vec_on_exit(true), myFirstID_(0), myNumIDs_(0), myCoefs_(NULL), nonlocalIDs_(NULL), nonlocalElementSize_(NULL), numNonlocalIDs_(0), allocatedNonlocalLength_(0), nonlocalCoefs_(NULL), last_edit(0), ignoreNonLocalEntries_(false) { this->_type = type; }
| libMesh::EpetraVector< T >::EpetraVector | ( | const Parallel::Communicator & | comm, |
| 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 622 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
: NumericVector<T>(comm, type), _destroy_vec_on_exit(true), myFirstID_(0), myNumIDs_(0), myCoefs_(NULL), nonlocalIDs_(NULL), nonlocalElementSize_(NULL), numNonlocalIDs_(0), allocatedNonlocalLength_(0), nonlocalCoefs_(NULL), last_edit(0), ignoreNonLocalEntries_(false) { this->init(n, n, false, type); }
| libMesh::EpetraVector< T >::EpetraVector | ( | const Parallel::Communicator & | comm, |
| 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 646 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
: NumericVector<T>(comm, type), _destroy_vec_on_exit(true), myFirstID_(0), myNumIDs_(0), myCoefs_(NULL), nonlocalIDs_(NULL), nonlocalElementSize_(NULL), numNonlocalIDs_(0), allocatedNonlocalLength_(0), nonlocalCoefs_(NULL), last_edit(0), ignoreNonLocalEntries_(false) { this->init(n, n_local, false, type); }
| libMesh::EpetraVector< T >::EpetraVector | ( | const Parallel::Communicator & | comm, |
| 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 711 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::init().
: NumericVector<T>(comm, AUTOMATIC), _destroy_vec_on_exit(true), myFirstID_(0), myNumIDs_(0), myCoefs_(NULL), nonlocalIDs_(NULL), nonlocalElementSize_(NULL), numNonlocalIDs_(0), allocatedNonlocalLength_(0), nonlocalCoefs_(NULL), last_edit(0), ignoreNonLocalEntries_(false) { this->init(n, n_local, ghost, false, type); }
| libMesh::EpetraVector< T >::EpetraVector | ( | Epetra_Vector & | v, |
| const Parallel::Communicator &comm | LIBMESH_CAN_DEFAULT_TO_COMMWORLD | ||
| ) |
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.
| libMesh::EpetraVector< T >::~EpetraVector | ( | ) | [inline] |
Destructor, deallocates memory. Made virtual to allow for derived classes to behave properly.
Definition at line 747 of file trilinos_epetra_vector.h.
{
this->clear ();
}
| void libMesh::EpetraVector< T >::abs | ( | ) | [virtual] |
v = abs(v)... that is, each entry in v is replaced by its absolute value.
Implements libMesh::NumericVector< T >.
Definition at line 301 of file trilinos_epetra_vector.C.
| void libMesh::EpetraVector< T >::add | ( | const numeric_index_type | i, |
| const T | value | ||
| ) | [virtual] |
v(i) += value
Implements libMesh::NumericVector< T >.
Definition at line 194 of file trilinos_epetra_vector.C.
{
int i = static_cast<int> (i_in);
T value = value_in;
libmesh_assert_less (i_in, this->size());
SumIntoGlobalValues(1, &i, &value);
this->_is_closed = false;
}
| void libMesh::EpetraVector< T >::add | ( | const T | s | ) | [virtual] |
. Addition of s to all components. Note that s is a scalar and not a vector.
Implements libMesh::NumericVector< T >.
Definition at line 249 of file trilinos_epetra_vector.C.
{
const unsigned int nl = _vec->MyLength();
T * values = _vec->Values();
for (unsigned int i=0; i<nl; i++)
values[i]+=v_in;
this->_is_closed = false;
}
| void libMesh::EpetraVector< T >::add | ( | const NumericVector< T > & | V | ) | [virtual] |
. Simple vector addition, equal to the operator +=.
Implements libMesh::NumericVector< T >.
Definition at line 263 of file trilinos_epetra_vector.C.
{
this->add (1., v);
}
| void libMesh::EpetraVector< T >::add | ( | const T | a, |
| const NumericVector< T > & | v | ||
| ) | [virtual] |
. Simple vector addition, equal to the operator +=.
Implements libMesh::NumericVector< T >.
Definition at line 270 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, and libMesh::EpetraVector< T >::size().
| void libMesh::EpetraVector< T >::add_vector | ( | const T * | v, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [virtual] |
where v is a pointer and each dof_indices[i] specifies where to add value v[i]
Reimplemented from libMesh::NumericVector< T >.
Definition at line 209 of file trilinos_epetra_vector.C.
{
libmesh_assert_equal_to (sizeof(numeric_index_type), sizeof(int));
SumIntoGlobalValues (dof_indices.size(),
(int*) &dof_indices[0],
const_cast<T*>(v));
}
| void libMesh::EpetraVector< T >::add_vector | ( | const NumericVector< T > & | V, |
| const SparseMatrix< T > & | A | ||
| ) | [virtual] |
, add the product of a SparseMatrix A and a NumericVector V to this, where this=U.
Implements libMesh::NumericVector< T >.
Definition at line 223 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, and libMesh::EpetraVector< T >::zero_clone().
{
const EpetraVector<T>* V = cast_ptr<const EpetraVector<T>*>(&V_in);
const EpetraMatrix<T>* A = cast_ptr<const EpetraMatrix<T>*>(&A_in);
// FIXME - does Trilinos let us do this *without* memory allocation?
UniquePtr<NumericVector<T> > temp = V->zero_clone();
EpetraVector<T>* tempV = cast_ptr<EpetraVector<T>*>(temp.get());
A->mat()->Multiply(false, *V->_vec, *tempV->_vec);
*this += *temp;
}
| void libMesh::NumericVector< T >::add_vector | ( | const std::vector< T > & | v, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [inline, inherited] |
where v is a std::vector and each dof_indices[i] specifies where to add value v[i]
Definition at line 806 of file numeric_vector.h.
References libMesh::libmesh_assert().
{
libmesh_assert(v.size() == dof_indices.size());
if (!v.empty())
this->add_vector(&v[0], dof_indices);
}
| void libMesh::NumericVector< T >::add_vector | ( | const NumericVector< T > & | V, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [virtual, inherited] |
where v is a NumericVector and each dof_indices[i] specifies where to add value v(i)
Definition at line 395 of file numeric_vector.C.
References libMesh::NumericVector< T >::size().
{
int n = dof_indices.size();
libmesh_assert_equal_to(v.size(), static_cast<unsigned>(n));
for (int i=0; i<n; i++)
this->add (dof_indices[i], v(i));
}
| void libMesh::NumericVector< T >::add_vector | ( | const DenseVector< T > & | V, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [inline, inherited] |
where v is a DenseVector and each dof_indices[i] specifies where to add value v(i)
Definition at line 818 of file numeric_vector.h.
References libMesh::DenseVector< T >::empty(), libMesh::libmesh_assert(), and libMesh::DenseVector< T >::size().
{
libmesh_assert(v.size() == dof_indices.size());
if (!v.empty())
this->add_vector(&v(0), dof_indices);
}
| void libMesh::NumericVector< T >::add_vector | ( | const NumericVector< T > & | v, |
| const ShellMatrix< T > & | a | ||
| ) | [inherited] |
, add the product of a ShellMatrix A and a NumericVector V to this, where this=U.
Definition at line 407 of file numeric_vector.C.
References libMesh::ShellMatrix< T >::vector_mult_add().
{
a.vector_mult_add(*this,v);
}
| void libMesh::EpetraVector< T >::add_vector_transpose | ( | const NumericVector< T > & | V, |
| const SparseMatrix< T > & | A_trans | ||
| ) | [virtual] |
, add the product of the transpose of a SparseMatrix A_trans and a NumericVector V to this, where this=U.
Implements libMesh::NumericVector< T >.
Definition at line 240 of file trilinos_epetra_vector.C.
{
libmesh_not_implemented();
}
| UniquePtr< NumericVector< T > > libMesh::NumericVector< T >::build | ( | const Parallel::Communicator & | comm, |
| const SolverPackage | solver_package = libMesh::default_solver_package() |
||
| ) | [static, inherited] |
Builds a NumericVector on the processors in communicator comm using the linear solver package specified by solver_package
Definition at line 46 of file numeric_vector.C.
References libMesh::AUTOMATIC, libMesh::EIGEN_SOLVERS, libMesh::LASPACK_SOLVERS, libMesh::PETSC_SOLVERS, and libMesh::TRILINOS_SOLVERS.
Referenced by libMesh::NumericVector< T >::build(), and libMesh::ExactErrorEstimator::estimate_error().
{
// Build the appropriate vector
switch (solver_package)
{
#ifdef LIBMESH_HAVE_LASPACK
case LASPACK_SOLVERS:
return UniquePtr<NumericVector<T> >(new LaspackVector<T>(comm, AUTOMATIC));
#endif
#ifdef LIBMESH_HAVE_PETSC
case PETSC_SOLVERS:
return UniquePtr<NumericVector<T> >(new PetscVector<T>(comm, AUTOMATIC));
#endif
#ifdef LIBMESH_HAVE_TRILINOS
case TRILINOS_SOLVERS:
return UniquePtr<NumericVector<T> >(new EpetraVector<T>(comm, AUTOMATIC));
#endif
#ifdef LIBMESH_HAVE_EIGEN
case EIGEN_SOLVERS:
return UniquePtr<NumericVector<T> >(new EigenSparseVector<T>(comm, AUTOMATIC));
#endif
default:
return UniquePtr<NumericVector<T> >(new DistributedVector<T>(comm, AUTOMATIC));
}
libmesh_error_msg("We'll never get here!");
return UniquePtr<NumericVector<T> >();
}
| UniquePtr< NumericVector< T > > libMesh::NumericVector< T >::build | ( | const SolverPackage | solver_package = libMesh::default_solver_package() | ) | [static, inherited] |
Builds a NumericVector on the processors in communicator CommWorld using the linear solver package specified by solver_package. Deprecated.
Definition at line 84 of file numeric_vector.C.
References libMesh::NumericVector< T >::build(), and libMesh::CommWorld.
{
libmesh_deprecated();
return NumericVector<T>::build(CommWorld, solver_package);
}
| void libMesh::EpetraVector< T >::clear | ( | ) | [inline, virtual] |
EpetraVector<T> to a pristine state. Reimplemented from libMesh::NumericVector< T >.
Definition at line 860 of file trilinos_epetra_vector.h.
References libMesh::libMeshPrivateData::_is_initialized, and libMesh::initialized().
{
if (this->initialized())
{
// We might just be an interface to a user-provided _vec
if (this->_destroy_vec_on_exit)
{
delete _vec;
_vec = NULL;
}
// But we currently always own our own _map
delete _map;
_map = NULL;
}
this->_is_closed = this->_is_initialized = false;
}
| UniquePtr< NumericVector< T > > libMesh::EpetraVector< T >::clone | ( | ) | const [inline, virtual] |
Creates a copy of this vector and returns it in an UniquePtr.
Implements libMesh::NumericVector< T >.
Definition at line 909 of file trilinos_epetra_vector.h.
References libMesh::AUTOMATIC.
| void libMesh::EpetraVector< T >::close | ( | ) | [inline, virtual] |
Call the assemble functions
Implements libMesh::NumericVector< T >.
Definition at line 836 of file trilinos_epetra_vector.h.
References libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::Parallel::Communicator::max().
Referenced by libMesh::EpetraVector< T >::localize().
{
libmesh_assert (this->initialized());
// Are we adding or inserting?
unsigned char global_last_edit = last_edit;
this->comm().max(global_last_edit);
libmesh_assert(!last_edit || last_edit == global_last_edit);
if (global_last_edit == 1)
this->GlobalAssemble(Insert);
else if (global_last_edit == 2)
this->GlobalAssemble(Add);
else
libmesh_assert(!global_last_edit);
this->_is_closed = true;
this->last_edit = 0;
}
| virtual bool libMesh::NumericVector< T >::closed | ( | ) | const [inline, virtual, inherited] |
Definition at line 147 of file numeric_vector.h.
Referenced by libMesh::DofMap::enforce_adjoint_constraints_exactly(), libMesh::DofMap::enforce_constraints_exactly(), libMesh::DofMap::max_constraint_error(), libMesh::LaspackVector< T >::operator=(), libMesh::EigenSparseVector< T >::operator=(), and libMesh::PetscVector< T >::operator=().
{ return _is_closed; }
| const Parallel::Communicator& libMesh::ParallelObject::comm | ( | ) | const [inline, inherited] |
Parallel::Communicator object used by this mesh. Definition at line 86 of file parallel_object.h.
References libMesh::ParallelObject::_communicator.
Referenced by libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_petsc_snes_residual(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::MetisPartitioner::_do_partition(), libMesh::ParmetisPartitioner::_do_repartition(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::bounding_box(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::MeshRefinement::coarsen_elements(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshSetSystem_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::ParmetisPartitioner::initialize(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::ParallelMesh::libmesh_assert_valid_parallel_flags(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::MeshSerializer::MeshSerializer(), libMesh::ParallelMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::MeshTools::n_p_levels(), libMesh::ParallelMesh::parallel_max_elem_id(), libMesh::ParallelMesh::parallel_max_node_id(), libMesh::ParallelMesh::parallel_n_elem(), libMesh::ParallelMesh::parallel_n_nodes(), libMesh::Partitioner::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::System::project_vector(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::MeshRefinement::refine_elements(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::Parallel::sync_element_data_by_parent_id(), libMesh::Parallel::sync_node_data_by_element_id(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::NameBasedIO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::LegacyXdrIO::write_mesh(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), and libMesh::DivaIO::write_stream().
{ return _communicator; }
| int libMesh::NumericVector< T >::compare | ( | const NumericVector< T > & | other_vector, |
| const Real | threshold = TOLERANCE |
||
| ) | const [virtual, inherited] |
-1 when this is equivalent to other_vector, up to the given threshold. When differences occur, the return value contains the first index i where the difference (a[i]-b[i]) exceeded the threshold. When no threshold is given, the libMesh TOLERANCE is used. Definition at line 116 of file numeric_vector.C.
References std::abs(), libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), and std::max().
{
libmesh_assert (this->initialized());
libmesh_assert (other_vector.initialized());
libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
int first_different_i = std::numeric_limits<int>::max();
numeric_index_type i = first_local_index();
do
{
if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
first_different_i = i;
else
i++;
}
while (first_different_i==std::numeric_limits<int>::max()
&& i<last_local_index());
// Find the correct first differing index in parallel
this->comm().min(first_different_i);
if (first_different_i == std::numeric_limits<int>::max())
return -1;
return first_different_i;
}
| void libMesh::EpetraVector< T >::conjugate | ( | ) | [virtual] |
Replace each entry v_i = real(v_i) + imag(v_i) of this vector by its complex conjugate, real(v_i) - imag(v_i). Epetra is real-valued only, rendering this a no-op.
Implements libMesh::NumericVector< T >.
Definition at line 186 of file trilinos_epetra_vector.C.
{
// EPetra is real, rendering this a no-op.
}
| 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.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 528 of file trilinos_epetra_vector.C.
{
libmesh_not_implemented();
}
| void libMesh::EpetraVector< T >::destroyNonlocalData | ( | ) | [private] |
Definition at line 854 of file trilinos_epetra_vector.C.
{
if (allocatedNonlocalLength_ > 0) {
delete [] nonlocalIDs_;
delete [] nonlocalElementSize_;
nonlocalIDs_ = NULL;
nonlocalElementSize_ = NULL;
for(int i=0; i<numNonlocalIDs_; ++i) {
delete [] nonlocalCoefs_[i];
}
delete [] nonlocalCoefs_;
nonlocalCoefs_ = NULL;
numNonlocalIDs_ = 0;
allocatedNonlocalLength_ = 0;
}
return;
}
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
Referenced by libMesh::LibMeshInit::LibMeshInit().
{
_enable_print_counter = false;
return;
}
| T libMesh::EpetraVector< T >::dot | ( | const NumericVector< T > & | V | ) | const [virtual] |
Computes the dot product, p = U.V
Implements libMesh::NumericVector< T >.
Definition at line 308 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
{
const EpetraVector<T>* V = cast_ptr<const EpetraVector<T>*>(&V_in);
T result=0.0;
_vec->Dot(*V->_vec, &result);
return result;
}
| virtual T libMesh::NumericVector< T >::el | ( | const numeric_index_type | i | ) | const [inline, virtual, inherited] |
| 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.
{
_enable_print_counter = true;
return;
}
| void libMesh::EpetraVector< T >::FEoperatorequals | ( | const EpetraVector< T > & | source | ) | [private] |
Definition at line 827 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::allocatedNonlocalLength_, libMesh::EpetraVector< T >::nonlocalCoefs_, libMesh::EpetraVector< T >::nonlocalElementSize_, libMesh::EpetraVector< T >::nonlocalIDs_, and libMesh::EpetraVector< T >::numNonlocalIDs_.
{
(*_vec) = *(source._vec);
destroyNonlocalData();
if (source.allocatedNonlocalLength_ > 0) {
allocatedNonlocalLength_ = source.allocatedNonlocalLength_;
numNonlocalIDs_ = source.numNonlocalIDs_;
nonlocalIDs_ = new int[allocatedNonlocalLength_];
nonlocalElementSize_ = new int[allocatedNonlocalLength_];
nonlocalCoefs_ = new double*[allocatedNonlocalLength_];
for(int i=0; i<numNonlocalIDs_; ++i) {
int elemSize = source.nonlocalElementSize_[i];
nonlocalCoefs_[i] = new double[elemSize];
nonlocalIDs_[i] = source.nonlocalIDs_[i];
nonlocalElementSize_[i] = elemSize;
for(int j=0; j<elemSize; ++j) {
nonlocalCoefs_[i][j] = source.nonlocalCoefs_[i][j];
}
}
}
}
| numeric_index_type libMesh::EpetraVector< T >::first_local_index | ( | ) | const [inline, virtual] |
Implements libMesh::NumericVector< T >.
Definition at line 945 of file trilinos_epetra_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
return _vec->Map().MinMyGID();
}
| void libMesh::NumericVector< T >::get | ( | const std::vector< numeric_index_type > & | index, |
| T * | values | ||
| ) | const [inline, virtual, inherited] |
Access multiple components at once. values will *not* be reallocated; it should already have enough space. The default implementation calls operator() for each index, but some implementations may supply faster methods here.
Reimplemented in libMesh::PetscVector< T >.
Definition at line 779 of file numeric_vector.h.
Referenced by libMesh::DofMap::enforce_adjoint_constraints_exactly(), libMesh::DofMap::enforce_constraints_exactly(), libMesh::FEMContext::pre_fe_reinit(), and libMesh::System::project_vector().
{
const std::size_t num = index.size();
for(std::size_t i=0; i<num; i++)
{
values[i] = (*this)(index[i]);
}
}
| void libMesh::NumericVector< T >::get | ( | const std::vector< numeric_index_type > & | index, |
| std::vector< T > & | values | ||
| ) | const [inline, 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.
Definition at line 792 of file numeric_vector.h.
{
const std::size_t num = index.size();
values.resize(num);
if (!num)
return;
this->get(index, &values[0]);
}
| 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().
{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
std::ostringstream oss;
oss << '\n'
<< " ---------------------------------------------------------------------------- \n"
<< "| Reference count information |\n"
<< " ---------------------------------------------------------------------------- \n";
for (Counts::iterator it = _counts.begin();
it != _counts.end(); ++it)
{
const std::string name(it->first);
const unsigned int creations = it->second.first;
const unsigned int destructions = it->second.second;
oss << "| " << name << " reference count information:\n"
<< "| Creations: " << creations << '\n'
<< "| Destructions: " << destructions << '\n';
}
oss << " ---------------------------------------------------------------------------- \n";
return oss.str();
#else
return "";
#endif
}
| int libMesh::NumericVector< T >::global_relative_compare | ( | const NumericVector< T > & | other_vector, |
| const Real | threshold = TOLERANCE |
||
| ) | const [virtual, inherited] |
-1 when this is equivalent to other_vector, up to the given local relative threshold. 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, the libMesh TOLERANCE is used. Definition at line 181 of file numeric_vector.C.
References std::abs(), libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), libMesh::NumericVector< T >::linfty_norm(), std::max(), and libMesh::Real.
{
libmesh_assert (this->initialized());
libmesh_assert (other_vector.initialized());
libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
int first_different_i = std::numeric_limits<int>::max();
numeric_index_type i = first_local_index();
const Real my_norm = this->linfty_norm();
const Real other_norm = other_vector.linfty_norm();
const Real abs_threshold = std::max(my_norm, other_norm) * threshold;
do
{
if ( std::abs( (*this)(i) - other_vector(i) ) > abs_threshold )
first_different_i = i;
else
i++;
}
while (first_different_i==std::numeric_limits<int>::max()
&& i<last_local_index());
// Find the correct first differing index in parallel
this->comm().min(first_different_i);
if (first_different_i == std::numeric_limits<int>::max())
return -1;
return first_different_i;
}
| int libMesh::EpetraVector< T >::GlobalAssemble | ( | Epetra_CombineMode | mode = Add | ) | [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 782 of file trilinos_epetra_vector.C.
{
//In this method we need to gather all the non-local (overlapping) data
//that's been input on each processor, into the (probably) non-overlapping
//distribution defined by the map that 'this' vector was constructed with.
//We don't need to do anything if there's only one processor or if
//ignoreNonLocalEntries_ is true.
if (_vec->Map().Comm().NumProc() < 2 || ignoreNonLocalEntries_) {
return(0);
}
//First build a map that describes the data in nonlocalIDs_/nonlocalCoefs_.
//We'll use the arbitrary distribution constructor of Map.
Epetra_BlockMap sourceMap(-1, numNonlocalIDs_,
nonlocalIDs_, nonlocalElementSize_,
_vec->Map().IndexBase(), _vec->Map().Comm());
//Now build a vector to hold our nonlocalCoefs_, and to act as the source-
//vector for our import operation.
Epetra_MultiVector nonlocalVector(sourceMap, 1);
int i,j;
for(i=0; i<numNonlocalIDs_; ++i) {
for(j=0; j<nonlocalElementSize_[i]; ++j) {
nonlocalVector.ReplaceGlobalValue(nonlocalIDs_[i], j, 0,
nonlocalCoefs_[i][j]);
}
}
Epetra_Export exporter(sourceMap, _vec->Map());
EPETRA_CHK_ERR( _vec->Export(nonlocalVector, exporter, mode) );
destroyNonlocalData();
return(0);
}
| 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, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.first++;
}
| 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, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.second++;
}
| 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 756 of file trilinos_epetra_vector.h.
References libMesh::libMeshPrivateData::_is_initialized, libMesh::AUTOMATIC, libMesh::GHOSTED, libMesh::libmesh_assert(), libMesh::PARALLEL, libMesh::SERIAL, and libMesh::zero.
Referenced by libMesh::EpetraVector< T >::EpetraVector(), and libMesh::EpetraVector< T >::localize().
{
// We default to allocating n_local local storage
numeric_index_type my_n_local = n_local;
if (type == AUTOMATIC)
{
if (n == n_local)
this->_type = SERIAL;
else
this->_type = PARALLEL;
}
else if (type == GHOSTED)
{
// We don't yet support GHOSTED Epetra vectors, so to get the
// same functionality we need a SERIAL vector with local
// storage allocated for every entry.
this->_type = SERIAL;
my_n_local = n;
}
else
this->_type = type;
libmesh_assert ((this->_type==SERIAL && n==my_n_local) ||
this->_type==PARALLEL);
_map = new Epetra_Map(static_cast<int>(n),
my_n_local,
0,
Epetra_MpiComm (this->comm().get()));
_vec = new Epetra_Vector(*_map);
myFirstID_ = _vec->Map().MinMyGID();
myNumIDs_ = _vec->Map().NumMyElements();
//Currently we impose the restriction that NumVectors==1, so we won't
//need the LDA argument when calling ExtractView. Hence the "dummy" arg.
int dummy;
_vec->ExtractView(&myCoefs_, &dummy);
this->_is_initialized = true;
this->_is_closed = true;
this->last_edit = 0;
if (fast == false)
this->zero ();
}
| 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 825 of file trilinos_epetra_vector.h.
References libMesh::TriangleWrapper::init().
| 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 811 of file trilinos_epetra_vector.h.
References libMesh::TriangleWrapper::init().
| void libMesh::EpetraVector< T >::init | ( | const NumericVector< T > & | other, |
| const bool | fast = false |
||
| ) | [virtual] |
Creates a vector that has the same dimension and storage type as other, including ghost dofs.
Implements libMesh::NumericVector< T >.
Definition at line 737 of file trilinos_epetra_vector.h.
References libMesh::TriangleWrapper::init(), libMesh::NumericVector< T >::local_size(), libMesh::NumericVector< T >::size(), and libMesh::NumericVector< T >::type().
{
this->init(other.size(),other.local_size(),fast,other.type());
}
| virtual bool libMesh::NumericVector< T >::initialized | ( | ) | const [inline, virtual, inherited] |
Definition at line 131 of file numeric_vector.h.
Referenced by libMesh::ImplicitSystem::assemble(), libMesh::NumericVector< T >::compare(), libMesh::PetscVector< T >::create_subvector(), libMesh::NumericVector< T >::global_relative_compare(), libMesh::PetscVector< T >::init(), and libMesh::NumericVector< T >::local_relative_compare().
{ return _is_initialized; }
| int libMesh::EpetraVector< T >::inputNonlocalValue | ( | int | GID, |
| double | value, | ||
| bool | accumulate | ||
| ) | [private] |
Definition at line 682 of file trilinos_epetra_vector.C.
{
int insertPoint = -1;
//find offset of GID in nonlocalIDs_
int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_,
insertPoint);
if (offset >= 0) {
//if offset >= 0
// put value in nonlocalCoefs_[offset][0]
if (accumulate) {
nonlocalCoefs_[offset][0] += value;
}
else {
nonlocalCoefs_[offset][0] = value;
}
}
else {
//else
// insert GID in nonlocalIDs_
// insert 1 in nonlocalElementSize_
// insert value in nonlocalCoefs_
int tmp1 = numNonlocalIDs_;
int tmp2 = allocatedNonlocalLength_;
int tmp3 = allocatedNonlocalLength_;
EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_,
tmp1, tmp2) );
--tmp1;
EPETRA_CHK_ERR( Epetra_Util_insert(1, insertPoint, nonlocalElementSize_,
tmp1, tmp3) );
double* values = new double[1];
values[0] = value;
EPETRA_CHK_ERR( Epetra_Util_insert(values, insertPoint, nonlocalCoefs_,
numNonlocalIDs_, allocatedNonlocalLength_) );
}
return(0);
}
| int libMesh::EpetraVector< T >::inputNonlocalValues | ( | int | GID, |
| int | numValues, | ||
| const double * | values, | ||
| bool | accumulate | ||
| ) | [private] |
Definition at line 725 of file trilinos_epetra_vector.C.
References libMesh::err.
{
int insertPoint = -1;
//find offset of GID in nonlocalIDs_
int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_,
insertPoint);
if (offset >= 0) {
//if offset >= 0
// put value in nonlocalCoefs_[offset][0]
if (numValues != nonlocalElementSize_[offset]) {
libMesh::err << "Epetra_FEVector ERROR: block-size for GID " << GID << " is "
<< numValues<<" which doesn't match previously set block-size of "
<< nonlocalElementSize_[offset] << std::endl;
return(-1);
}
if (accumulate) {
for(int j=0; j<numValues; ++j) {
nonlocalCoefs_[offset][j] += values[j];
}
}
else {
for(int j=0; j<numValues; ++j) {
nonlocalCoefs_[offset][j] = values[j];
}
}
}
else {
//else
// insert GID in nonlocalIDs_
// insert numValues in nonlocalElementSize_
// insert values in nonlocalCoefs_
int tmp1 = numNonlocalIDs_;
int tmp2 = allocatedNonlocalLength_;
int tmp3 = allocatedNonlocalLength_;
EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_,
tmp1, tmp2) );
--tmp1;
EPETRA_CHK_ERR( Epetra_Util_insert(numValues, insertPoint, nonlocalElementSize_,
tmp1, tmp3) );
double* newvalues = new double[numValues];
for(int j=0; j<numValues; ++j) {
newvalues[j] = values[j];
}
EPETRA_CHK_ERR( Epetra_Util_insert(newvalues, insertPoint, nonlocalCoefs_,
numNonlocalIDs_, allocatedNonlocalLength_) );
}
return(0);
}
| int libMesh::EpetraVector< T >::inputValues | ( | int | numIDs, |
| const int * | GIDs, | ||
| const double * | values, | ||
| bool | accumulate | ||
| ) | [private] |
Definition at line 602 of file trilinos_epetra_vector.C.
References libMesh::libmesh_assert().
{
if (accumulate) {
libmesh_assert(last_edit == 0 || last_edit == 2);
last_edit = 2;
} else {
libmesh_assert(last_edit == 0 || last_edit == 1);
last_edit = 1;
}
//Important note!! This method assumes that there is only 1 point
//associated with each element.
for(int i=0; i<numIDs; ++i) {
if (_vec->Map().MyGID(GIDs[i])) {
if (accumulate) {
_vec->SumIntoGlobalValue(GIDs[i], 0, 0, values[i]);
}
else {
_vec->ReplaceGlobalValue(GIDs[i], 0, 0, values[i]);
}
}
else {
if (!ignoreNonLocalEntries_) {
EPETRA_CHK_ERR( inputNonlocalValue(GIDs[i], values[i], accumulate) );
}
}
}
return(0);
}
| int libMesh::EpetraVector< T >::inputValues | ( | int | numIDs, |
| const int * | GIDs, | ||
| const int * | numValuesPerID, | ||
| const double * | values, | ||
| bool | accumulate | ||
| ) | [private] |
Definition at line 639 of file trilinos_epetra_vector.C.
References libMesh::libmesh_assert().
{
if (accumulate) {
libmesh_assert(last_edit == 0 || last_edit == 2);
last_edit = 2;
} else {
libmesh_assert(last_edit == 0 || last_edit == 1);
last_edit = 1;
}
int offset=0;
for(int i=0; i<numIDs; ++i) {
int numValues = numValuesPerID[i];
if (_vec->Map().MyGID(GIDs[i])) {
if (accumulate) {
for(int j=0; j<numValues; ++j) {
_vec->SumIntoGlobalValue(GIDs[i], j, 0, values[offset+j]);
}
}
else {
for(int j=0; j<numValues; ++j) {
_vec->ReplaceGlobalValue(GIDs[i], j, 0, values[offset+j]);
}
}
}
else {
if (!ignoreNonLocalEntries_) {
EPETRA_CHK_ERR( inputNonlocalValues(GIDs[i], numValues,
&(values[offset]), accumulate) );
}
}
offset += numValues;
}
return(0);
}
| void libMesh::EpetraVector< T >::insert | ( | const T * | v, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [virtual] |
where v is a T[] or T* and you want to specify WHERE to insert it
Reimplemented from libMesh::NumericVector< T >.
Definition at line 282 of file trilinos_epetra_vector.C.
{
libmesh_assert_equal_to (sizeof(numeric_index_type), sizeof(int));
ReplaceGlobalValues (dof_indices.size(),
(int*) &dof_indices[0],
const_cast<T*>(v));
}
| void libMesh::NumericVector< T >::insert | ( | const std::vector< T > & | v, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [inline, inherited] |
where v is a std::vector<T> and you want to specify WHERE to insert it
Definition at line 830 of file numeric_vector.h.
References libMesh::libmesh_assert().
{
libmesh_assert(v.size() == dof_indices.size());
if (!v.empty())
this->insert(&v[0], dof_indices);
}
| void libMesh::NumericVector< T >::insert | ( | const NumericVector< T > & | V, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [virtual, inherited] |
, where U and V are type NumericVector<T> and you want to specify WHERE to insert the NumericVector<T> V
Definition at line 104 of file numeric_vector.C.
References libMesh::NumericVector< T >::size().
{
libmesh_assert_equal_to (V.size(), dof_indices.size());
for (numeric_index_type i=0; i<dof_indices.size(); i++)
this->set (dof_indices[i], V(i));
}
| void libMesh::NumericVector< T >::insert | ( | const DenseVector< T > & | V, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [inline, inherited] |
where U and V are type DenseVector<T> and you want to specify WHERE to insert the DenseVector<T> V
Definition at line 842 of file numeric_vector.h.
References libMesh::DenseVector< T >::empty(), libMesh::libmesh_assert(), and libMesh::DenseVector< T >::size().
{
libmesh_assert(v.size() == dof_indices.size());
if (!v.empty())
this->insert(&v(0), dof_indices);
}
| void libMesh::NumericVector< T >::insert | ( | const DenseSubVector< T > & | V, |
| const std::vector< numeric_index_type > & | dof_indices | ||
| ) | [inline, inherited] |
where V is a DenseSubVector<T> and you want to specify WHERE to insert it
Definition at line 854 of file numeric_vector.h.
References libMesh::DenseVectorBase< T >::empty(), libMesh::libmesh_assert(), and libMesh::DenseVectorBase< T >::size().
{
libmesh_assert(v.size() == dof_indices.size());
if (!v.empty())
this->insert(&v(0), dof_indices);
}
| Real libMesh::EpetraVector< T >::l1_norm | ( | ) | const [virtual] |
-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::closed(), libMesh::libmesh_assert(), and libMesh::Real.
{
libmesh_assert(this->closed());
Real value;
_vec->Norm1(&value);
return value;
}
| Real libMesh::EpetraVector< T >::l2_norm | ( | ) | const [virtual] |
-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::closed(), libMesh::libmesh_assert(), and libMesh::Real.
{
libmesh_assert(this->closed());
Real value;
_vec->Norm2(&value);
return value;
}
| numeric_index_type libMesh::EpetraVector< T >::last_local_index | ( | ) | const [inline, virtual] |
Implements libMesh::NumericVector< T >.
Definition at line 956 of file trilinos_epetra_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
return _vec->Map().MaxMyGID()+1;
}
| Real libMesh::EpetraVector< T >::linfty_norm | ( | ) | const [virtual] |
-norm of a vector. Implements libMesh::NumericVector< T >.
Definition at line 91 of file trilinos_epetra_vector.C.
References libMesh::closed(), libMesh::libmesh_assert(), and libMesh::Real.
{
libmesh_assert(this->closed());
Real value;
_vec->NormInf(&value);
return value;
}
| int libMesh::NumericVector< T >::local_relative_compare | ( | const NumericVector< T > & | other_vector, |
| const Real | threshold = TOLERANCE |
||
| ) | const [virtual, inherited] |
-1 when this is equivalent to other_vector, up to the given local relative threshold. 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, the libMesh TOLERANCE is used. Definition at line 148 of file numeric_vector.C.
References std::abs(), libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), and std::max().
{
libmesh_assert (this->initialized());
libmesh_assert (other_vector.initialized());
libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
int first_different_i = std::numeric_limits<int>::max();
numeric_index_type i = first_local_index();
do
{
if ( std::abs( (*this)(i) - other_vector(i) ) > threshold *
std::max(std::abs((*this)(i)), std::abs(other_vector(i))))
first_different_i = i;
else
i++;
}
while (first_different_i==std::numeric_limits<int>::max()
&& i<last_local_index());
// Find the correct first differing index in parallel
this->comm().min(first_different_i);
if (first_different_i == std::numeric_limits<int>::max())
return -1;
return first_different_i;
}
| numeric_index_type libMesh::EpetraVector< T >::local_size | ( | ) | const [inline, virtual] |
Implements libMesh::NumericVector< T >.
Definition at line 936 of file trilinos_epetra_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
return _vec->MyLength();
}
| void libMesh::EpetraVector< T >::localize | ( | std::vector< T > & | v_local | ) | const [virtual] |
Creates a copy of the global vector in the local vector v_local.
Implements libMesh::NumericVector< T >.
Definition at line 479 of file trilinos_epetra_vector.C.
References libMesh::libmesh_assert().
Referenced by libMesh::EpetraVector< T >::localize().
{
// This function must be run on all processors at once
parallel_object_only();
const unsigned int n = this->size();
const unsigned int nl = this->local_size();
libmesh_assert(this->_vec);
v_local.clear();
v_local.reserve(n);
// build up my local part
for (unsigned int i=0; i<nl; i++)
v_local.push_back((*this->_vec)[i]);
this->comm().allgather (v_local);
}
| void libMesh::EpetraVector< T >::localize | ( | NumericVector< T > & | v_local | ) | const [virtual] |
Same, but fills a NumericVector<T> instead of a std::vector.
Implements libMesh::NumericVector< T >.
Definition at line 407 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
| void libMesh::EpetraVector< T >::localize | ( | NumericVector< T > & | v_local, |
| const std::vector< numeric_index_type > & | send_list | ||
| ) | const [virtual] |
Creates a local vector v_local containing only information relevant to this processor, as defined by the send_list.
Implements libMesh::NumericVector< T >.
Definition at line 421 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::localize().
Referenced by libMesh::EpetraVector< T >::localize().
{
// TODO: optimize to sync only the send list values
this->localize(v_local_in);
// EpetraVector<T>* v_local =
// cast_ptr<EpetraVector<T>*>(&v_local_in);
// libmesh_assert(this->_map.get());
// libmesh_assert(v_local->_map.get());
// libmesh_assert_equal_to (v_local->local_size(), this->size());
// libmesh_assert_less_equal (send_list.size(), v_local->size());
// Epetra_Import importer (*v_local->_map, *this->_map);
// v_local->_vec->Import (*this->_vec, importer, Insert);
}
| 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 >.
Definition at line 442 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::close(), libMesh::EpetraVector< T >::init(), libMesh::EpetraVector< T >::localize(), libMesh::PARALLEL, and libMesh::EpetraVector< T >::set().
{
// Only good for serial vectors.
libmesh_assert_equal_to (this->size(), this->local_size());
libmesh_assert_greater (last_local_idx, first_local_idx);
libmesh_assert_less_equal (send_list.size(), this->size());
libmesh_assert_less (last_local_idx, this->size());
const unsigned int my_size = this->size();
const unsigned int my_local_size = (last_local_idx - first_local_idx + 1);
// Don't bother for serial cases
if ((first_local_idx == 0) &&
(my_local_size == my_size))
return;
// Build a parallel vector, initialize it with the local
// parts of (*this)
EpetraVector<T> parallel_vec(this->comm(), PARALLEL);
parallel_vec.init (my_size, my_local_size, true, PARALLEL);
// Copy part of *this into the parallel_vec
for (numeric_index_type i=first_local_idx; i<=last_local_idx; i++)
parallel_vec.set(i,this->el(i));
// localize like normal
parallel_vec.close();
parallel_vec.localize (*this, send_list);
this->close();
}
| void libMesh::EpetraVector< T >::localize_to_one | ( | std::vector< T > & | v_local, |
| const processor_id_type | proc_id = 0 |
||
| ) | const [virtual] |
Creates a local copy of the global vector in v_local only on processor proc_id. By default the data is sent to processor 0. This method is useful for outputting data from one processor.
Implements libMesh::NumericVector< T >.
Definition at line 502 of file trilinos_epetra_vector.C.
References libMesh::libmesh_assert(), and libMesh::n_processors().
{
// This function must be run on all processors at once
parallel_object_only();
const unsigned int n = this->size();
const unsigned int nl = this->local_size();
libmesh_assert_less (pid, this->n_processors());
libmesh_assert(this->_vec);
v_local.clear();
v_local.reserve(n);
// build up my local part
for (unsigned int i=0; i<nl; i++)
v_local.push_back((*this->_vec)[i]);
this->comm().gather (pid, v_local);
}
| Real libMesh::EpetraVector< T >::max | ( | ) | const [inline, virtual] |
Implements libMesh::NumericVector< T >.
Definition at line 994 of file trilinos_epetra_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
T value;
_vec->MaxValue(&value);
return value;
}
| Real libMesh::EpetraVector< T >::min | ( | ) | const [inline, virtual] |
Implements libMesh::NumericVector< T >.
Definition at line 979 of file trilinos_epetra_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
T value;
_vec->MinValue(&value);
return value;
}
| 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.
Referenced by libMesh::LibMeshInit::~LibMeshInit().
{ return _n_objects; }
| processor_id_type libMesh::ParallelObject::n_processors | ( | ) | const [inline, inherited] |
Definition at line 92 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::size().
Referenced by libMesh::ParmetisPartitioner::_do_repartition(), libMesh::ParallelMesh::add_elem(), libMesh::ParallelMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::ParallelMesh::assign_unique_ids(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::ParallelMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshBase::get_info(), libMesh::EquationSystems::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshTools::processor_bounding_box(), libMesh::System::project_vector(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::Partitioner::repartition(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::BoundaryInfo::sync(), libMesh::ParallelMesh::update_parallel_id_counts(), libMesh::CheckpointIO::write(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().
{ return cast_int<processor_id_type>(_communicator.size()); }
| 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 966 of file trilinos_epetra_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
libmesh_assert ( ((i >= this->first_local_index()) &&
(i < this->last_local_index())) );
return (*_vec)[i-this->first_local_index()];
}
| NumericVector<T>& libMesh::NumericVector< T >::operator*= | ( | const T | a | ) | [inline, inherited] |
Multiplication operator. Equivalent to U.scale(a)
Definition at line 378 of file numeric_vector.h.
{ this->scale(a); return *this; }
| NumericVector< T > & libMesh::EpetraVector< T >::operator+= | ( | const NumericVector< T > & | V | ) | [virtual] |
Addition operator. Fast equivalent to U.add(1, V).
Implements libMesh::NumericVector< T >.
Definition at line 104 of file trilinos_epetra_vector.C.
References libMesh::closed(), and libMesh::libmesh_assert().
{
libmesh_assert(this->closed());
this->add(1., v);
return *this;
}
| NumericVector< T > & libMesh::EpetraVector< T >::operator-= | ( | const NumericVector< T > & | V | ) | [virtual] |
Subtraction operator. Fast equivalent to U.add(-1, V).
Implements libMesh::NumericVector< T >.
Definition at line 117 of file trilinos_epetra_vector.C.
References libMesh::closed(), and libMesh::libmesh_assert().
{
libmesh_assert(this->closed());
this->add(-1., v);
return *this;
}
| NumericVector< T > & libMesh::EpetraVector< T >::operator/= | ( | NumericVector< T > & | v | ) | [virtual] |
Pointwise Division operator. ie divide every entry in this vector by the entry in v
Implements libMesh::NumericVector< T >.
Definition at line 129 of file trilinos_epetra_vector.C.
References libMesh::closed(), libMesh::libmesh_assert(), and libMesh::NumericVector< T >::size().
{
libmesh_assert(this->closed());
libmesh_assert_equal_to(size(), v.size());
EpetraVector<T> & v_vec = cast_ref<EpetraVector<T>&>(v);
_vec->ReciprocalMultiply(1.0, *v_vec._vec, *_vec, 0.0);
}
| NumericVector<T>& libMesh::NumericVector< T >::operator/= | ( | const T | a | ) | [inline, inherited] |
Division operator. Equivalent to U.scale(1./a)
Definition at line 384 of file numeric_vector.h.
{ this->scale(1./a); return *this; }
| NumericVector< T > & libMesh::EpetraVector< T >::operator= | ( | const T | s | ) | [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 333 of file trilinos_epetra_vector.C.
{
_vec->PutScalar(s_in);
return *this;
}
| NumericVector< T > & libMesh::EpetraVector< T >::operator= | ( | const NumericVector< T > & | V | ) | [virtual] |
: copy all components.
Implements libMesh::NumericVector< T >.
Definition at line 344 of file trilinos_epetra_vector.C.
{
const EpetraVector<T>* v = cast_ptr<const EpetraVector<T>*>(&v_in);
*this = *v;
return *this;
}
| EpetraVector< T > & libMesh::EpetraVector< T >::operator= | ( | const EpetraVector< T > & | V | ) |
: copy all components.
Definition at line 357 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
{
(*_vec) = *v._vec;
// FIXME - what about our communications data?
return *this;
}
| NumericVector< T > & libMesh::EpetraVector< T >::operator= | ( | const std::vector< T > & | v | ) | [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 370 of file trilinos_epetra_vector.C.
{
T * values = _vec->Values();
if(this->size() == v.size())
{
const unsigned int nl=this->local_size();
const unsigned int fli=this->first_local_index();
for(unsigned int i=0;i<nl;i++)
values[i]=v[fli+i];
}
else
{
libmesh_assert_equal_to (v.size(), this->local_size());
const unsigned int nl=this->local_size();
for(unsigned int i=0;i<nl;i++)
values[i]=v[i];
}
return *this;
}
| void libMesh::EpetraVector< T >::pointwise_mult | ( | const NumericVector< T > & | vec1, |
| const NumericVector< T > & | vec2 | ||
| ) | [virtual] |
Computes the pointwise (i.e. component-wise) product of vec1 and vec2 and stores the result in *this.
Implements libMesh::NumericVector< T >.
Definition at line 321 of file trilinos_epetra_vector.C.
References libMesh::EpetraVector< T >::_vec.
{
const EpetraVector<T>* V1 = cast_ptr<const EpetraVector<T>*>(&vec1);
const EpetraVector<T>* V2 = cast_ptr<const EpetraVector<T>*>(&vec2);
_vec->Multiply(1.0, *V1->_vec, *V2->_vec, 0.0);
}
| 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 887 of file numeric_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
os << "Size\tglobal = " << this->size()
<< "\t\tlocal = " << this->local_size() << std::endl;
os << "#\tValue" << std::endl;
for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++)
os << i << "\t" << (*this)(i) << std::endl;
}
| void libMesh::NumericVector< Complex >::print | ( | std::ostream & | os | ) | const [inline, inherited] |
Definition at line 869 of file numeric_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
os << "Size\tglobal = " << this->size()
<< "\t\tlocal = " << this->local_size() << std::endl;
// std::complex<>::operator<<() is defined, but use this form
os << "#\tReal part\t\tImaginary part" << std::endl;
for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++)
os << i << "\t"
<< (*this)(i).real() << "\t\t"
<< (*this)(i).imag() << std::endl;
}
| 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 924 of file numeric_vector.h.
References libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::processor_id().
{
libmesh_assert (this->initialized());
std::vector<T> v(this->size());
this->localize(v);
// Right now we only want one copy of the output
if (this->processor_id())
return;
os << "Size\tglobal = " << this->size() << std::endl;
os << "#\tValue" << std::endl;
for (numeric_index_type i=0; i!=v.size(); i++)
os << i << "\t" << v[i] << std::endl;
}
| void libMesh::NumericVector< Complex >::print_global | ( | std::ostream & | os | ) | const [inline, inherited] |
Definition at line 902 of file numeric_vector.h.
References libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::processor_id().
{
libmesh_assert (this->initialized());
std::vector<Complex> v(this->size());
this->localize(v);
// Right now we only want one copy of the output
if (this->processor_id())
return;
os << "Size\tglobal = " << this->size() << std::endl;
os << "#\tReal part\t\tImaginary part" << std::endl;
for (numeric_index_type i=0; i!=v.size(); i++)
os << i << "\t"
<< v[i].real() << "\t\t"
<< v[i].imag() << std::endl;
}
| 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().
Referenced by libMesh::LibMeshInit::~LibMeshInit().
{
if( _enable_print_counter ) out_stream << ReferenceCounter::get_info();
}
| virtual void libMesh::NumericVector< T >::print_matlab | ( | const std::string & | = "" | ) | const [inline, virtual, inherited] |
Print the contents of the matrix in Matlab's sparse matrix format. Optionally prints the matrix to the file named name. If name is not specified it is dumped to the screen.
Reimplemented in libMesh::PetscVector< T >.
Definition at line 645 of file numeric_vector.h.
{
libmesh_not_implemented();
}
| processor_id_type libMesh::ParallelObject::processor_id | ( | ) | const [inline, inherited] |
Definition at line 98 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::rank().
Referenced by libMesh::MetisPartitioner::_do_partition(), libMesh::EquationSystems::_read_impl(), libMesh::SerialMesh::active_local_elements_begin(), libMesh::ParallelMesh::active_local_elements_begin(), libMesh::SerialMesh::active_local_elements_end(), libMesh::ParallelMesh::active_local_elements_end(), libMesh::SerialMesh::active_local_subdomain_elements_begin(), libMesh::ParallelMesh::active_local_subdomain_elements_begin(), libMesh::SerialMesh::active_local_subdomain_elements_end(), libMesh::ParallelMesh::active_local_subdomain_elements_end(), libMesh::SerialMesh::active_not_local_elements_begin(), libMesh::ParallelMesh::active_not_local_elements_begin(), libMesh::SerialMesh::active_not_local_elements_end(), libMesh::ParallelMesh::active_not_local_elements_end(), libMesh::ParallelMesh::add_elem(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::ParallelMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::ParallelMesh::assign_unique_ids(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::ParmetisPartitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::DofMap::build_sparsity(), libMesh::ParallelMesh::clear(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::ExodusII_IO_Helper::create(), libMesh::ParallelMesh::delete_elem(), libMesh::ParallelMesh::delete_node(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::SerialMesh::facelocal_elements_begin(), libMesh::ParallelMesh::facelocal_elements_begin(), libMesh::SerialMesh::facelocal_elements_end(), libMesh::ParallelMesh::facelocal_elements_end(), libMesh::MeshFunction::find_element(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::ParallelMesh::insert_elem(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::SerialMesh::local_elements_begin(), libMesh::ParallelMesh::local_elements_begin(), libMesh::SerialMesh::local_elements_end(), libMesh::ParallelMesh::local_elements_end(), libMesh::SerialMesh::local_level_elements_begin(), libMesh::ParallelMesh::local_level_elements_begin(), libMesh::SerialMesh::local_level_elements_end(), libMesh::ParallelMesh::local_level_elements_end(), libMesh::SerialMesh::local_nodes_begin(), libMesh::ParallelMesh::local_nodes_begin(), libMesh::SerialMesh::local_nodes_end(), libMesh::ParallelMesh::local_nodes_end(), libMesh::SerialMesh::local_not_level_elements_begin(), libMesh::ParallelMesh::local_not_level_elements_begin(), libMesh::SerialMesh::local_not_level_elements_end(), libMesh::ParallelMesh::local_not_level_elements_end(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::SerialMesh::not_local_elements_begin(), libMesh::ParallelMesh::not_local_elements_begin(), libMesh::SerialMesh::not_local_elements_end(), libMesh::ParallelMesh::not_local_elements_end(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::ParallelMesh::ParallelMesh(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::System::project_vector(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::MeshData::read_xdr(), libMesh::SerialMesh::semilocal_elements_begin(), libMesh::ParallelMesh::semilocal_elements_begin(), libMesh::SerialMesh::semilocal_elements_end(), libMesh::ParallelMesh::semilocal_elements_end(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::BoundaryInfo::sync(), libMesh::MeshTools::total_weight(), libMesh::ParallelMesh::update_parallel_id_counts(), libMesh::MeshTools::weight(), libMesh::NameBasedIO::write(), libMesh::ExodusII_IO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), and libMesh::ExodusII_IO_Helper::write_timestep().
{ return cast_int<processor_id_type>(_communicator.rank()); }
| void libMesh::EpetraVector< T >::reciprocal | ( | ) | [virtual] |
Replace each entry v_i of this vector by its reciprocal, 1/v_i.
Implements libMesh::NumericVector< T >.
Definition at line 158 of file trilinos_epetra_vector.C.
References std::abs(), and std::min().
{
// The Epetra::reciprocal() function takes a constant reference to *another* vector,
// and fills _vec with its reciprocal. Does that mean we can't pass *_vec as the
// argument?
// _vec->reciprocal( *_vec );
// Alternatively, compute the reciprocal by hand... see also the add(T) member that does this...
const unsigned int nl = _vec->MyLength();
T* values = _vec->Values();
for (unsigned int i=0; i<nl; i++)
{
// Don't divide by zero (maybe only check this in debug mode?)
if (std::abs(values[i]) < std::numeric_limits<T>::min())
libmesh_error_msg("Error, divide by zero in DistributedVector<T>::reciprocal()!");
values[i] = 1. / values[i];
}
// Leave the vector in a closed state...
this->close();
}
| int libMesh::EpetraVector< T >::ReplaceGlobalValues | ( | int | numIDs, |
| const int * | GIDs, | ||
| const double * | values | ||
| ) | [private] |
Copy values into the vector overwriting any values that already exist for the specified indices.
Definition at line 573 of file trilinos_epetra_vector.C.
{
return( inputValues( numIDs, GIDs, values, false) );
}
| int libMesh::EpetraVector< T >::ReplaceGlobalValues | ( | const Epetra_IntSerialDenseVector & | GIDs, |
| const Epetra_SerialDenseVector & | values | ||
| ) | [private] |
Copy values into the vector, replacing any values that already exist for the specified GIDs.
| 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 581 of file trilinos_epetra_vector.C.
{
if (GIDs.Length() != values.Length()) {
return(-1);
}
return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), false) );
}
| int libMesh::EpetraVector< T >::ReplaceGlobalValues | ( | int | numIDs, |
| const int * | GIDs, | ||
| const int * | numValuesPerID, | ||
| const double * | values | ||
| ) | [private] |
Definition at line 593 of file trilinos_epetra_vector.C.
{
return( inputValues( numIDs, GIDs, numValuesPerID, values, false) );
}
| void libMesh::EpetraVector< T >::scale | ( | const T | factor | ) | [virtual] |
Scale each element of the vector by the given factor.
Implements libMesh::NumericVector< T >.
Definition at line 295 of file trilinos_epetra_vector.C.
{
_vec->Scale(factor_in);
}
| void libMesh::EpetraVector< T >::set | ( | const numeric_index_type | i, |
| const T | value | ||
| ) | [virtual] |
v(i) = value
Implements libMesh::NumericVector< T >.
Definition at line 143 of file trilinos_epetra_vector.C.
Referenced by libMesh::EpetraVector< T >::localize().
{
int i = static_cast<int> (i_in);
T value = value_in;
libmesh_assert_less (i_in, this->size());
ReplaceGlobalValues(1, &i, &value);
this->_is_closed = false;
}
| void libMesh::EpetraVector< T >::setIgnoreNonLocalEntries | ( | bool | flag | ) | [inline, private] |
Set whether or not non-local data values should be ignored.
Definition at line 551 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::ignoreNonLocalEntries_.
{
ignoreNonLocalEntries_ = flag;
}
| numeric_index_type libMesh::EpetraVector< T >::size | ( | ) | const [inline, virtual] |
n(), but was renamed to get the EpetraVector<T> class closer to the C++ standard library's std::vector container. Implements libMesh::NumericVector< T >.
Definition at line 925 of file trilinos_epetra_vector.h.
References libMesh::initialized(), and libMesh::libmesh_assert().
Referenced by libMesh::EpetraVector< T >::add().
{
libmesh_assert (this->initialized());
return _vec->GlobalLength();
}
| Real libMesh::NumericVector< T >::subset_l1_norm | ( | const std::set< numeric_index_type > & | indices | ) | const [virtual, inherited] |
-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 324 of file numeric_vector.C.
References std::abs(), and libMesh::Real.
Referenced by libMesh::System::discrete_var_norm().
| Real libMesh::NumericVector< T >::subset_l2_norm | ( | const std::set< numeric_index_type > & | indices | ) | const [virtual, inherited] |
-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 342 of file numeric_vector.C.
References libMesh::TensorTools::norm_sq(), and libMesh::Real.
Referenced by libMesh::System::discrete_var_norm().
{
const NumericVector<T> & v = *this;
std::set<numeric_index_type>::const_iterator it = indices.begin();
const std::set<numeric_index_type>::const_iterator it_end = indices.end();
Real norm = 0;
for(; it!=it_end; ++it)
norm += TensorTools::norm_sq(v(*it));
this->comm().sum(norm);
return std::sqrt(norm);
}
| Real libMesh::NumericVector< T >::subset_linfty_norm | ( | const std::set< numeric_index_type > & | indices | ) | const [virtual, inherited] |
-norm of a vector.Note that the indices must necessary live on this processor.
Definition at line 360 of file numeric_vector.C.
References std::abs(), and libMesh::Real.
Referenced by libMesh::System::discrete_var_norm().
{
const NumericVector<T> & v = *this;
std::set<numeric_index_type>::const_iterator it = indices.begin();
const std::set<numeric_index_type>::const_iterator it_end = indices.end();
Real norm = 0;
for(; it!=it_end; ++it)
{
Real value = std::abs(v(*it));
if(value > norm)
norm = value;
}
this->comm().max(norm);
return norm;
}
| T libMesh::EpetraVector< T >::sum | ( | ) | const [virtual] |
Implements libMesh::NumericVector< T >.
Definition at line 48 of file trilinos_epetra_vector.C.
References libMesh::closed(), libMesh::libmesh_assert(), and libMesh::Parallel::sum().
| int libMesh::EpetraVector< T >::SumIntoGlobalValues | ( | int | numIDs, |
| const int * | GIDs, | ||
| const double * | values | ||
| ) | [private] |
Accumulate values into the vector, adding them to any values that already exist for the specified indices.
Definition at line 544 of file trilinos_epetra_vector.C.
{
return( inputValues( numIDs, GIDs, values, true) );
}
| int libMesh::EpetraVector< T >::SumIntoGlobalValues | ( | const Epetra_IntSerialDenseVector & | GIDs, |
| const Epetra_SerialDenseVector & | values | ||
| ) | [private] |
Accumulate values into the vector, adding them to any values that already exist for the specified GIDs.
| 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 552 of file trilinos_epetra_vector.C.
{
if (GIDs.Length() != values.Length()) {
return(-1);
}
return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), true) );
}
| int libMesh::EpetraVector< T >::SumIntoGlobalValues | ( | int | numIDs, |
| const int * | GIDs, | ||
| const int * | numValuesPerID, | ||
| const double * | values | ||
| ) | [private] |
Definition at line 564 of file trilinos_epetra_vector.C.
{
return( inputValues( numIDs, GIDs, numValuesPerID, values, true) );
}
| void libMesh::EpetraVector< T >::swap | ( | NumericVector< T > & | v | ) | [inline, virtual] |
Swaps the raw Epetra vector context pointers.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 1009 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_destroy_vec_on_exit, libMesh::EpetraVector< T >::_map, libMesh::EpetraVector< T >::_vec, libMesh::EpetraVector< T >::allocatedNonlocalLength_, libMesh::EpetraVector< T >::ignoreNonLocalEntries_, libMesh::EpetraVector< T >::last_edit, libMesh::EpetraVector< T >::myCoefs_, libMesh::EpetraVector< T >::myFirstID_, libMesh::EpetraVector< T >::myNumIDs_, libMesh::EpetraVector< T >::nonlocalCoefs_, libMesh::EpetraVector< T >::nonlocalElementSize_, libMesh::EpetraVector< T >::nonlocalIDs_, libMesh::EpetraVector< T >::numNonlocalIDs_, libMesh::swap(), and libMesh::NumericVector< T >::swap().
{
NumericVector<T>::swap(other);
EpetraVector<T>& v = cast_ref<EpetraVector<T>&>(other);
std::swap(_vec, v._vec);
std::swap(_map, v._map);
std::swap(_destroy_vec_on_exit, v._destroy_vec_on_exit);
std::swap(myFirstID_, v.myFirstID_);
std::swap(myNumIDs_, v.myNumIDs_);
std::swap(myCoefs_, v.myCoefs_);
std::swap(nonlocalIDs_, v.nonlocalIDs_);
std::swap(nonlocalElementSize_, v.nonlocalElementSize_);
std::swap(numNonlocalIDs_, v.numNonlocalIDs_);
std::swap(allocatedNonlocalLength_, v.allocatedNonlocalLength_);
std::swap(nonlocalCoefs_, v.nonlocalCoefs_);
std::swap(last_edit, v.last_edit);
std::swap(ignoreNonLocalEntries_, v.ignoreNonLocalEntries_);
}
| ParallelType libMesh::NumericVector< T >::type | ( | ) | const [inline, inherited] |
Definition at line 136 of file numeric_vector.h.
Referenced by libMesh::DofMap::enforce_adjoint_constraints_exactly(), libMesh::DofMap::enforce_constraints_exactly(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::MeshFunction::find_element(), libMesh::DistributedVector< T >::init(), libMesh::LaspackVector< T >::init(), libMesh::EigenSparseVector< T >::init(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< T >::localize(), libMesh::PetscVector< T >::operator=(), libMesh::System::project_vector(), and libMesh::System::read_serialized_vector().
{ return _type; }
| ParallelType& libMesh::NumericVector< T >::type | ( | ) | [inline, inherited] |
Definition at line 141 of file numeric_vector.h.
{ 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 468 of file trilinos_epetra_vector.h.
References libMesh::EpetraVector< T >::_vec, and libMesh::libmesh_assert().
Referenced by libMesh::EpetraMatrix< T >::get_diagonal(), and libMesh::NoxNonlinearSolver< T >::solve().
{ libmesh_assert(_vec); return _vec; }
| 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 883 of file trilinos_epetra_vector.h.
References libMesh::closed(), libMesh::initialized(), and libMesh::libmesh_assert().
{
libmesh_assert (this->initialized());
libmesh_assert (this->closed());
_vec->PutScalar(0.0);
}
| UniquePtr< 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 UniquePtr.
Implements libMesh::NumericVector< T >.
Definition at line 895 of file trilinos_epetra_vector.h.
References libMesh::AUTOMATIC.
Referenced by libMesh::EpetraVector< T >::add_vector().
| 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 633 of file numeric_vector.h.
{
v.print_global(os);
return os;
}
const Parallel::Communicator& libMesh::ParallelObject::_communicator [protected, inherited] |
Definition at line 104 of file parallel_object.h.
Referenced by libMesh::EquationSystems::build_solution_vector(), libMesh::ParallelObject::comm(), libMesh::EquationSystems::get_solution(), libMesh::ParallelObject::n_processors(), libMesh::ParallelObject::operator=(), and libMesh::ParallelObject::processor_id().
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 487 of file trilinos_epetra_vector.h.
Referenced by 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 675 of file numeric_vector.h.
Referenced by libMesh::NumericVector< Number >::closed(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::operator=(), libMesh::EigenSparseVector< T >::swap(), 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 681 of file numeric_vector.h.
Referenced by libMesh::PetscVector< T >::create_subvector(), libMesh::NumericVector< Number >::initialized(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::operator=(), libMesh::EigenSparseVector< T >::swap(), and libMesh::NumericVector< T >::swap().
Epetra_Map* libMesh::EpetraVector< T >::_map [private] |
Holds the distributed Map
Definition at line 481 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::swap().
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 686 of file numeric_vector.h.
Referenced by libMesh::DistributedVector< T >::DistributedVector(), libMesh::EigenSparseVector< T >::EigenSparseVector(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::PetscVector< T >::init(), libMesh::LaspackVector< T >::LaspackVector(), libMesh::PetscVector< T >::operator=(), libMesh::PetscVector< T >::PetscVector(), libMesh::EigenSparseVector< T >::swap(), 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 476 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::add(), libMesh::EpetraVector< T >::add_vector(), libMesh::EpetraVector< T >::dot(), libMesh::EpetraVector< T >::FEoperatorequals(), libMesh::EpetraVector< T >::localize(), libMesh::EpetraVector< T >::operator=(), libMesh::EpetraVector< T >::pointwise_mult(), libMesh::EpetraVector< T >::swap(), and libMesh::EpetraVector< T >::vec().
int libMesh::EpetraVector< T >::allocatedNonlocalLength_ [private] |
Definition at line 580 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::FEoperatorequals(), and libMesh::EpetraVector< T >::swap().
bool libMesh::EpetraVector< T >::ignoreNonLocalEntries_ [private] |
Definition at line 590 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::setIgnoreNonLocalEntries(), and libMesh::EpetraVector< T >::swap().
unsigned char libMesh::EpetraVector< T >::last_edit [private] |
Keep track of whether the last write operation on this vector was nothing (0) or a sum (1) or an add (2), so we can decide how to do the GlobalAssemble()
Definition at line 588 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::swap().
double* libMesh::EpetraVector< T >::myCoefs_ [private] |
Definition at line 575 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::swap().
int libMesh::EpetraVector< T >::myFirstID_ [private] |
Definition at line 573 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::swap().
int libMesh::EpetraVector< T >::myNumIDs_ [private] |
Definition at line 574 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::swap().
double** libMesh::EpetraVector< T >::nonlocalCoefs_ [private] |
Definition at line 581 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::FEoperatorequals(), and libMesh::EpetraVector< T >::swap().
int* libMesh::EpetraVector< T >::nonlocalElementSize_ [private] |
Definition at line 578 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::FEoperatorequals(), and libMesh::EpetraVector< T >::swap().
int* libMesh::EpetraVector< T >::nonlocalIDs_ [private] |
Definition at line 577 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::FEoperatorequals(), and libMesh::EpetraVector< T >::swap().
int libMesh::EpetraVector< T >::numNonlocalIDs_ [private] |
Definition at line 579 of file trilinos_epetra_vector.h.
Referenced by libMesh::EpetraVector< T >::FEoperatorequals(), and libMesh::EpetraVector< T >::swap().