A vector in the mathematical sense, initially developed for linear algebra type applications. If MPI then this vector can be distributed - its distribution is described by the LinearAlgebraDistribution object at Distribution_pt. Data is stored in a C-style pointer vector (double*) More...
#include <double_vector.h>
Public Member Functions | |
DoubleVector () | |
Constructor for an uninitialized DoubleVector. More... | |
DoubleVector (const LinearAlgebraDistribution *const &dist_pt, const double &v=0.0) | |
Constructor. Assembles a DoubleVector with a prescribed distribution. Additionally every entry can be set (with argument v - defaults to 0). More... | |
DoubleVector (const LinearAlgebraDistribution &dist, const double &v=0.0) | |
Constructor. Assembles a DoubleVector with a prescribed distribution. Additionally every entry can be set (with argument v - defaults to 0). More... | |
~DoubleVector () | |
Destructor - just calls this->clear() to delete the distribution and data. More... | |
DoubleVector (const DoubleVector &new_vector) | |
Copy constructor. More... | |
void | operator= (const DoubleVector &old_vector) |
assignment operator More... | |
void | build (const DoubleVector &old_vector) |
Just copys the argument DoubleVector. More... | |
void | build (const LinearAlgebraDistribution &dist, const double &v) |
Assembles a DoubleVector with distribution dist, if v is specified each element is set to v, otherwise each element is set to 0.0. More... | |
void | build (const LinearAlgebraDistribution *const &dist_pt, const double &v) |
Assembles a DoubleVector with distribution dist, if v is specified each element is set to v, otherwise each element is set to 0.0. More... | |
void | build (const LinearAlgebraDistribution &dist, const Vector< double > &v) |
Assembles a DoubleVector with a distribution dist and coefficients taken from the vector v. Note. The vector v MUST be of length nrow() More... | |
void | build (const LinearAlgebraDistribution *const &dist_pt, const Vector< double > &v) |
Assembles a DoubleVector with a distribution dist and coefficients taken from the vector v. Note. The vector v MUST be of length nrow() More... | |
void | initialise (const double &v) |
initialise the whole vector with value v More... | |
void | initialise (const Vector< double > v) |
initialise the vector with coefficient from the vector v. Note: The vector v must be of length More... | |
void | clear () |
wipes the DoubleVector More... | |
bool | built () const |
void | set_external_values (const LinearAlgebraDistribution *const &dist_pt, double *external_values, bool delete_external_values) |
Allows are external data to be used by this vector. WARNING: The size of the external data must correspond to the LinearAlgebraDistribution dist_pt argument. More... | |
void | set_external_values (double *external_values, bool delete_external_values) |
Allows are external data to be used by this vector. WARNING: The size of the external data must correspond to the distribution of this vector. More... | |
void | redistribute (const LinearAlgebraDistribution *const &dist_pt) |
The contents of the vector are redistributed to match the new distribution. In a non-MPI rebuild this method works, but does nothing. NOTE 1: The current distribution and the new distribution must have the same number of global rows. NOTE 2: The current distribution and the new distribution must have the same Communicator. More... | |
double & | operator[] (int i) |
[] access function to the (local) values of this vector More... | |
bool | operator== (const DoubleVector &v) |
== operator More... | |
void | operator+= (const DoubleVector &v) |
+= operator with another vector More... | |
void | operator-= (const DoubleVector &v) |
-= operator with another vector More... | |
void | operator*= (const double &d) |
multiply by a double More... | |
void | operator/= (const double &d) |
divide by a double More... | |
const double & | operator[] (int i) const |
[] access function to the (local) values of this vector More... | |
double | max () const |
returns the maximum coefficient More... | |
double * | values_pt () |
access function to the underlying values More... | |
double * | values_pt () const |
access function to the underlying values (const version) More... | |
void | output (std::ostream &outfile, const int &output_precision=-1) const |
output the global contents of the vector More... | |
void | output (std::string filename, const int &output_precision=-1) const |
output the global contents of the vector More... | |
void | output_local_values (std::ostream &outfile, const int &output_precision=-1) const |
output the local contents of the vector More... | |
void | output_local_values (std::string filename, const int &output_precision=-1) const |
output the local contents of the vector More... | |
void | output_local_values_with_offset (std::ostream &outfile, const int &output_precision=-1) const |
output the local contents of the vector More... | |
void | output_local_values_with_offset (std::string filename, const int &output_precision=-1) const |
output the local contents of the vector More... | |
double | dot (const DoubleVector &vec) const |
compute the dot product of this vector with the vector vec. More... | |
double | norm () const |
compute the 2 norm of this vector More... | |
double | norm (const CRDoubleMatrix *matrix_pt) const |
compute the A-norm using the matrix at matrix_pt More... | |
Public Member Functions inherited from oomph::DistributableLinearAlgebraObject | |
DistributableLinearAlgebraObject () | |
Default constructor - create a distribution. More... | |
DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix)=delete | |
Broken copy constructor. More... | |
void | operator= (const DistributableLinearAlgebraObject &)=delete |
Broken assignment operator. More... | |
virtual | ~DistributableLinearAlgebraObject () |
Destructor. More... | |
LinearAlgebraDistribution * | distribution_pt () const |
access to the LinearAlgebraDistribution More... | |
unsigned | nrow () const |
access function to the number of global rows. More... | |
unsigned | nrow_local () const |
access function for the num of local rows on this processor. More... | |
unsigned | nrow_local (const unsigned &p) const |
access function for the num of local rows on this processor. More... | |
unsigned | first_row () const |
access function for the first row on this processor More... | |
unsigned | first_row (const unsigned &p) const |
access function for the first row on this processor More... | |
bool | distributed () const |
distribution is serial or distributed More... | |
bool | distribution_built () const |
if the communicator_pt is null then the distribution is not setup then false is returned, otherwise return true More... | |
void | build_distribution (const LinearAlgebraDistribution *const dist_pt) |
setup the distribution of this distributable linear algebra object More... | |
void | build_distribution (const LinearAlgebraDistribution &dist) |
setup the distribution of this distributable linear algebra object More... | |
Private Attributes | |
double * | Values_pt |
the local vector More... | |
bool | Internal_values |
Boolean flag to indicate whether the vector's data (values_pt) is owned by this vector. More... | |
bool | Built |
indicates that the vector has been built and is usable More... | |
Friends | |
std::ostream & | operator<< (std::ostream &out, const DoubleVector &v) |
Ouput operator for DoubleVector. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from oomph::DistributableLinearAlgebraObject | |
void | clear_distribution () |
clear the distribution of this distributable linear algebra object More... | |
A vector in the mathematical sense, initially developed for linear algebra type applications. If MPI then this vector can be distributed - its distribution is described by the LinearAlgebraDistribution object at Distribution_pt. Data is stored in a C-style pointer vector (double*)
Definition at line 57 of file double_vector.h.
|
inline |
Constructor for an uninitialized DoubleVector.
Definition at line 61 of file double_vector.h.
|
inline |
Constructor. Assembles a DoubleVector with a prescribed distribution. Additionally every entry can be set (with argument v - defaults to 0).
Definition at line 66 of file double_vector.h.
References build().
|
inline |
Constructor. Assembles a DoubleVector with a prescribed distribution. Additionally every entry can be set (with argument v - defaults to 0).
Definition at line 76 of file double_vector.h.
References build().
|
inline |
Destructor - just calls this->clear() to delete the distribution and data.
Definition at line 84 of file double_vector.h.
References clear().
|
inline |
void oomph::DoubleVector::build | ( | const DoubleVector & | old_vector | ) |
Just copys the argument DoubleVector.
Definition at line 35 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::distribution_built(), oomph::DistributableLinearAlgebraObject::distribution_pt(), Internal_values, oomph::DistributableLinearAlgebraObject::nrow_local(), values_pt(), and Values_pt.
Referenced by oomph::OomphLibPreconditionerEpetraOperator::ApplyInverse(), build(), oomph::ComplexDampedJacobi< MATRIX >::complex_solve_helper(), oomph::DoubleVectorHelpers::concatenate(), oomph::DoubleVectorHelpers::concatenate_without_communication(), DoubleVector(), oomph::BlockPreconditioner< MATRIX >::get_block_vector(), oomph::Problem::get_dofs(), oomph::Problem::get_inverse_mass_matrix_times_residuals(), oomph::Problem::get_jacobian(), oomph::Problem::get_residuals(), oomph::HypreInterface::hypre_solve(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_ordered_preconditioner_vector(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_vector(), oomph::TrilinosEpetraHelpers::multiply(), oomph::CRDoubleMatrix::multiply(), oomph::DenseDoubleMatrix::multiply(), oomph::CCDoubleMatrix::multiply(), oomph::MatrixVectorProduct::multiply(), oomph::CRDoubleMatrix::multiply_transpose(), oomph::DenseDoubleMatrix::multiply_transpose(), oomph::CCDoubleMatrix::multiply_transpose(), oomph::MatrixVectorProduct::multiply_transpose(), operator=(), oomph::DoubleVectorWithHaloEntries::operator=(), oomph::PitchForkHandler::PitchForkHandler(), oomph::BlockDiagonalPreconditioner< MATRIX >::preconditioner_solve(), oomph::DummyBlockPreconditioner< MATRIX >::preconditioner_solve(), oomph::MatrixBasedDiagPreconditioner::preconditioner_solve(), oomph::MatrixBasedLumpedPreconditioner< MATRIX >::preconditioner_solve(), oomph::TrilinosPreconditionerBase::preconditioner_solve(), oomph::FSIPreconditioner::preconditioner_solve(), oomph::PseudoElasticFSIPreconditioner::preconditioner_solve(), oomph::LagrangeEnforcedFlowPreconditioner::preconditioner_solve(), oomph::NavierStokesSchurComplementPreconditioner::preconditioner_solve(), oomph::PressureBasedSolidLSCPreconditioner::preconditioner_solve(), oomph::SpaceTimeNavierStokesSubsidiaryPreconditioner::preconditioner_solve(), oomph::AugmentedBlockFoldLinearSolver::resolve(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::AugmentedBlockPitchForkLinearSolver::resolve(), oomph::HSL_MA42::resolve(), oomph::TrilinosAztecOOSolver::resolve(), oomph::Problem::setup_element_count_per_dof(), oomph::DenseLU::solve(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::HSL_MA42::solve(), oomph::CG< MATRIX >::solve(), oomph::BiCGStab< MATRIX >::solve(), oomph::GS< CRDoubleMatrix >::solve(), oomph::GMRES< MATRIX >::solve(), oomph::AugmentedProblemGMRES::solve(), oomph::SuperLUSolver::solve(), oomph::MumpsSolver::solve(), oomph::HelmholtzGMRESMG< MATRIX >::solve(), oomph::HelmholtzFGMRESMG< MATRIX >::solve(), oomph::TrilinosAztecOOSolver::solve(), oomph::BlockHopfLinearSolver::solve_for_two_rhs(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), oomph::GS< MATRIX >::solve_helper(), oomph::GS< CRDoubleMatrix >::solve_helper(), oomph::DampedJacobi< MATRIX >::solve_helper(), oomph::GMRES< MATRIX >::solve_helper(), oomph::AugmentedProblemGMRES::solve_helper(), oomph::HelmholtzGMRESMG< MATRIX >::solve_helper(), oomph::HelmholtzFGMRESMG< MATRIX >::solve_helper(), oomph::SuperLUSolver::solve_transpose(), and oomph::LowStorageRungeKutta< ORDER >::timestep().
|
inline |
Assembles a DoubleVector with distribution dist, if v is specified each element is set to v, otherwise each element is set to 0.0.
Definition at line 110 of file double_vector.h.
References build().
|
inline |
Assembles a DoubleVector with a distribution dist and coefficients taken from the vector v. Note. The vector v MUST be of length nrow()
Definition at line 123 of file double_vector.h.
References build().
void oomph::DoubleVector::build | ( | const LinearAlgebraDistribution *const & | dist_pt, |
const double & | v | ||
) |
Assembles a DoubleVector with distribution dist, if v is specified each element is set to v, otherwise each element is set to 0.0.
Assembles a DoubleVector with distribution dist, if v is specified each row is set to v.
Definition at line 59 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::build_distribution(), Built, oomph::LinearAlgebraDistribution::built(), clear(), Internal_values, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
void oomph::DoubleVector::build | ( | const LinearAlgebraDistribution *const & | dist_pt, |
const Vector< double > & | v | ||
) |
Assembles a DoubleVector with a distribution dist and coefficients taken from the vector v. Note. The vector v MUST be of length nrow()
Definition at line 91 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::build_distribution(), Built, oomph::LinearAlgebraDistribution::built(), clear(), initialise(), Internal_values, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
|
inline |
Definition at line 154 of file double_vector.h.
References Built.
Referenced by oomph::MumpsSolver::backsub(), oomph::SuperLUSolver::backsub_serial(), oomph::SuperLUSolver::backsub_transpose_serial(), oomph::Smoother::check_validity_of_solve_helper_inputs(), oomph::DoubleVectorHelpers::concatenate(), oomph::DoubleVectorHelpers::concatenate_without_communication(), oomph::TrilinosEpetraHelpers::copy_to_oomphlib_vector(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_vector(), oomph::TrilinosEpetraHelpers::create_epetra_vector_view_data(), dot(), oomph::BlockPreconditioner< MATRIX >::get_block_ordered_preconditioner_vector(), oomph::BlockPreconditioner< MATRIX >::get_block_vector(), oomph::BlockPreconditioner< MATRIX >::get_block_vectors(), oomph::Problem::get_jacobian(), oomph::Problem::get_residuals(), oomph::HypreInterface::hypre_solve(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_ordered_preconditioner_vector(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_vector(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_vectors(), oomph::BlockPreconditioner< MATRIX >::internal_return_block_ordered_preconditioner_vector(), oomph::BlockPreconditioner< MATRIX >::internal_return_block_vector(), oomph::BlockPreconditioner< MATRIX >::internal_return_block_vectors(), oomph::CRDoubleMatrix::lubksub(), oomph::TrilinosEpetraHelpers::multiply(), oomph::CRDoubleMatrix::multiply(), oomph::DenseDoubleMatrix::multiply(), oomph::CCDoubleMatrix::multiply(), oomph::MatrixVectorProduct::multiply(), oomph::CRDoubleMatrix::multiply_transpose(), oomph::DenseDoubleMatrix::multiply_transpose(), oomph::CCDoubleMatrix::multiply_transpose(), oomph::MatrixVectorProduct::multiply_transpose(), norm(), operator+=(), operator-=(), operator==(), oomph::BlockDiagonalPreconditioner< MATRIX >::preconditioner_solve(), oomph::MatrixBasedDiagPreconditioner::preconditioner_solve(), oomph::MatrixBasedLumpedPreconditioner< MATRIX >::preconditioner_solve(), oomph::ILUZeroPreconditioner< CCDoubleMatrix >::preconditioner_solve(), oomph::ILUZeroPreconditioner< CRDoubleMatrix >::preconditioner_solve(), oomph::HyprePreconditioner::preconditioner_solve(), oomph::IdentityPreconditioner::preconditioner_solve(), oomph::FSIPreconditioner::preconditioner_solve(), oomph::LagrangeEnforcedFlowPreconditioner::preconditioner_solve(), oomph::NavierStokesSchurComplementPreconditioner::preconditioner_solve(), oomph::PressureBasedSolidLSCPreconditioner::preconditioner_solve(), oomph::AugmentedBlockFoldLinearSolver::resolve(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::AugmentedBlockPitchForkLinearSolver::resolve(), oomph::HypreSolver::resolve(), oomph::TrilinosAztecOOSolver::resolve(), oomph::BlockPreconditioner< MATRIX >::return_block_ordered_preconditioner_vector(), oomph::BlockPreconditioner< MATRIX >::return_block_vector(), oomph::BlockPreconditioner< MATRIX >::return_block_vectors(), oomph::BlockPreconditioner< MATRIX >::return_concatenated_block_vector(), oomph::MumpsSolver::solve(), oomph::HypreSolver::solve(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::BiCGStab< MATRIX >::solve(), oomph::GMRES< MATRIX >::solve(), oomph::AugmentedProblemGMRES::solve(), oomph::FD_LU::solve(), oomph::SuperLUSolver::solve(), oomph::HelmholtzGMRESMG< MATRIX >::solve(), oomph::HelmholtzFGMRESMG< MATRIX >::solve(), oomph::TrilinosAztecOOSolver::solve(), oomph::BlockHopfLinearSolver::solve_for_two_rhs(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), oomph::GMRES< MATRIX >::solve_helper(), oomph::AugmentedProblemGMRES::solve_helper(), oomph::HelmholtzGMRESMG< MATRIX >::solve_helper(), oomph::HelmholtzFGMRESMG< MATRIX >::solve_helper(), oomph::SuperLUSolver::solve_transpose(), oomph::DoubleVectorHelpers::split(), and oomph::DoubleVectorHelpers::split_without_communication().
|
inline |
wipes the DoubleVector
Definition at line 142 of file double_vector.h.
References Built, oomph::DistributableLinearAlgebraObject::clear_distribution(), Internal_values, and Values_pt.
Referenced by build(), oomph::Problem::newton_solve(), oomph::Problem::newton_solve_continuation(), oomph::BiharmonicPreconditioner::preconditioner_solve(), oomph::InexactSubBiharmonicPreconditioner::preconditioner_solve(), oomph::FSIPreconditioner::preconditioner_solve(), oomph::PseudoElasticFSIPreconditioner::preconditioner_solve(), oomph::LagrangeEnforcedFlowPreconditioner::preconditioner_solve(), oomph::NavierStokesSchurComplementPreconditioner::preconditioner_solve(), oomph::PressureBasedSolidLSCPreconditioner::preconditioner_solve(), oomph::SpaceTimeNavierStokesSubsidiaryPreconditioner::preconditioner_solve(), oomph::PseudoElasticPreconditionerSubsidiaryBlockPreconditionerOld::preconditioner_solve(), oomph::LinearSolver::reset_gradient(), set_external_values(), oomph::ARPACK::solve_eigenproblem_legacy(), and ~DoubleVector().
double oomph::DoubleVector::dot | ( | const DoubleVector & | vec | ) | const |
compute the dot product of this vector with the vector vec.
compute the dot product of this vector with the vector vec
Definition at line 805 of file double_vector.cc.
References built(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DistributableLinearAlgebraObject::nrow_local(), values_pt(), and Values_pt.
Referenced by oomph::AugmentedProblemGMRES::augmented_matrix_multiply(), oomph::Problem::calculate_continuation_derivatives_helper(), norm(), oomph::GMRESBlockPreconditioner::preconditioner_solve(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), and oomph::AugmentedProblemGMRES::solve_helper().
void oomph::DoubleVector::initialise | ( | const double & | v | ) |
initialise the whole vector with value v
Definition at line 123 of file double_vector.cc.
References Built, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by build(), oomph::Problem::get_jacobian(), oomph::CRDoubleMatrix::multiply(), oomph::DenseDoubleMatrix::multiply(), oomph::CCDoubleMatrix::multiply(), oomph::CRDoubleMatrix::multiply_transpose(), oomph::DenseDoubleMatrix::multiply_transpose(), oomph::CCDoubleMatrix::multiply_transpose(), oomph::BiharmonicPreconditioner::preconditioner_solve(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), oomph::GS< MATRIX >::solve_helper(), oomph::GS< CRDoubleMatrix >::solve_helper(), oomph::DampedJacobi< MATRIX >::solve_helper(), oomph::GMRES< MATRIX >::solve_helper(), oomph::AugmentedProblemGMRES::solve_helper(), oomph::HelmholtzGMRESMG< MATRIX >::solve_helper(), and oomph::HelmholtzFGMRESMG< MATRIX >::solve_helper().
void oomph::DoubleVector::initialise | ( | const Vector< double > | v | ) |
initialise the vector with coefficient from the vector v. Note: The vector v must be of length
Definition at line 138 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::first_row(), oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
double oomph::DoubleVector::max | ( | ) | const |
returns the maximum coefficient
Definition at line 604 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DistributableLinearAlgebraObject::nrow(), oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by oomph::DoubleMatrixBase::max_residual(), oomph::Problem::newton_solve(), oomph::Problem::newton_solve_continuation(), oomph::SegregatableFSIProblem::segregated_solve(), and oomph::SolidICProblem::set_newmark_initial_condition_consistently().
double oomph::DoubleVector::norm | ( | ) | const |
compute the 2 norm of this vector
Definition at line 867 of file double_vector.cc.
References built(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by oomph::ComplexDampedJacobi< MATRIX >::complex_solve_helper(), oomph::HypreInterface::hypre_solve(), oomph::PitchForkHandler::PitchForkHandler(), oomph::GMRESBlockPreconditioner::preconditioner_solve(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), oomph::GS< MATRIX >::solve_helper(), oomph::GS< CRDoubleMatrix >::solve_helper(), oomph::DampedJacobi< MATRIX >::solve_helper(), oomph::GMRES< MATRIX >::solve_helper(), and oomph::AugmentedProblemGMRES::solve_helper().
double oomph::DoubleVector::norm | ( | const CRDoubleMatrix * | matrix_pt | ) | const |
compute the A-norm using the matrix at matrix_pt
Definition at line 914 of file double_vector.cc.
References built(), oomph::CRDoubleMatrix::built(), oomph::DistributableLinearAlgebraObject::distribution_pt(), dot(), and oomph::CRDoubleMatrix::multiply().
void oomph::DoubleVector::operator*= | ( | const double & | d | ) |
multiply by a double
Multiply by double.
Definition at line 550 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::distribution_built(), i, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by operator/=().
void oomph::DoubleVector::operator+= | ( | const DoubleVector & | v | ) |
+= operator with another vector
+= operator
Definition at line 459 of file double_vector.cc.
References built(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DistributableLinearAlgebraObject::nrow_local(), values_pt(), and Values_pt.
void oomph::DoubleVector::operator-= | ( | const DoubleVector & | v | ) |
-= operator with another vector
-= operator
Definition at line 504 of file double_vector.cc.
References built(), oomph::DistributableLinearAlgebraObject::distribution_built(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DistributableLinearAlgebraObject::nrow_local(), values_pt(), and Values_pt.
void oomph::DoubleVector::operator/= | ( | const double & | d | ) |
divide by a double
Divide by double.
Definition at line 573 of file double_vector.cc.
References operator*=().
|
inline |
bool oomph::DoubleVector::operator== | ( | const DoubleVector & | v | ) |
== operator
Definition at line 426 of file double_vector.cc.
References built(), i, oomph::DistributableLinearAlgebraObject::nrow_local(), values_pt(), and Values_pt.
double & oomph::DoubleVector::operator[] | ( | int | i | ) |
[] access function to the (local) values of this vector
Definition at line 408 of file double_vector.cc.
References i, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
const double & oomph::DoubleVector::operator[] | ( | int | i | ) | const |
[] access function to the (local) values of this vector
Definition at line 586 of file double_vector.cc.
References i, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
void oomph::DoubleVector::output | ( | std::ostream & | outfile, |
const int & | output_precision = -1 |
||
) | const |
output the global contents of the vector
output the contents of the vector
Definition at line 653 of file double_vector.cc.
References oomph::LinearAlgebraDistribution::communicator_pt(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_pt(), oomph::DistributableLinearAlgebraObject::first_row(), i, oomph::DistributableLinearAlgebraObject::nrow(), oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by output().
|
inline |
output the global contents of the vector
Definition at line 269 of file double_vector.h.
References output().
void oomph::DoubleVector::output_local_values | ( | std::ostream & | outfile, |
const int & | output_precision = -1 |
||
) | const |
output the local contents of the vector
Definition at line 742 of file double_vector.cc.
References i, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by output_local_values().
|
inline |
output the local contents of the vector
Definition at line 283 of file double_vector.h.
References output_local_values().
void oomph::DoubleVector::output_local_values_with_offset | ( | std::ostream & | outfile, |
const int & | output_precision = -1 |
||
) | const |
output the local contents of the vector
output the local contents of the vector with the first row offset.
Definition at line 772 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::first_row(), i, oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by output_local_values_with_offset().
|
inline |
output the local contents of the vector
Definition at line 298 of file double_vector.h.
References output_local_values_with_offset().
void oomph::DoubleVector::redistribute | ( | const LinearAlgebraDistribution *const & | dist_pt | ) |
The contents of the vector are redistributed to match the new distribution. In a non-MPI rebuild this method works, but does nothing. NOTE 1: The current distribution and the new distribution must have the same number of global rows. NOTE 2: The current distribution and the new distribution must have the same Communicator.
The contents of the vector are redistributed to match the new distribution. In a non-MPI build this method works, but does nothing. NOTE 1: The current distribution and the new distribution must have the same number of global rows. NOTE 2: The current distribution and the new distribution must have the same Communicator.
Definition at line 164 of file double_vector.cc.
References oomph::DistributableLinearAlgebraObject::build_distribution(), oomph::LinearAlgebraDistribution::communicator_pt(), oomph::LinearAlgebraDistribution::distributed(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_pt(), oomph::LinearAlgebraDistribution::first_row(), oomph::DistributableLinearAlgebraObject::first_row(), i, Internal_values, oomph::LinearAlgebraDistribution::nrow(), oomph::DistributableLinearAlgebraObject::nrow(), oomph::LinearAlgebraDistribution::nrow_local(), oomph::DistributableLinearAlgebraObject::nrow_local(), and Values_pt.
Referenced by oomph::Problem::add_eigenvector_to_dofs(), oomph::Problem::assign_eigenvector_to_dofs(), oomph::MumpsSolver::backsub(), oomph::Problem::calculate_continuation_derivatives_helper(), oomph::Problem::get_jacobian(), oomph::HypreInterface::hypre_solve(), oomph::Problem::newton_solve(), oomph::Problem::newton_solve_continuation(), oomph::ILUZeroPreconditioner< CCDoubleMatrix >::preconditioner_solve(), oomph::ILUZeroPreconditioner< CRDoubleMatrix >::preconditioner_solve(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::CG< MATRIX >::solve(), oomph::BiCGStab< MATRIX >::solve(), oomph::GMRES< MATRIX >::solve(), oomph::AugmentedProblemGMRES::solve(), oomph::SuperLUSolver::solve(), oomph::MumpsSolver::solve(), oomph::HelmholtzGMRESMG< MATRIX >::solve(), oomph::HelmholtzFGMRESMG< MATRIX >::solve(), oomph::TrilinosAztecOOSolver::solve(), and oomph::SuperLUSolver::solve_transpose().
|
inline |
Allows are external data to be used by this vector. WARNING: The size of the external data must correspond to the LinearAlgebraDistribution dist_pt argument.
Definition at line 167 of file double_vector.h.
References oomph::DistributableLinearAlgebraObject::build_distribution(), Built, and clear().
Referenced by oomph::OomphLibPreconditionerEpetraOperator::ApplyInverse(), oomph::Problem::get_jacobian(), and oomph::Problem::get_residuals().
|
inline |
Allows are external data to be used by this vector. WARNING: The size of the external data must correspond to the distribution of this vector.
Definition at line 191 of file double_vector.h.
References oomph::DistributableLinearAlgebraObject::distribution_built(), Internal_values, and Values_pt.
|
inline |
access function to the underlying values
Definition at line 254 of file double_vector.h.
References Values_pt.
Referenced by oomph::DenseLU::backsub(), oomph::SuperLUSolver::backsub_distributed(), oomph::SuperLUSolver::backsub_serial(), oomph::SuperLUSolver::backsub_transpose_serial(), build(), oomph::Problem::calculate_continuation_derivatives_helper(), oomph::DoubleVectorHelpers::concatenate(), oomph::DoubleVectorHelpers::concatenate_without_communication(), oomph::TrilinosEpetraHelpers::copy_to_oomphlib_vector(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_vector(), oomph::TrilinosEpetraHelpers::create_epetra_vector_view_data(), oomph::HypreHelpers::create_HYPRE_Vector(), dot(), oomph::HypreInterface::hypre_solve(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_ordered_preconditioner_vector(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_vector(), oomph::BlockPreconditioner< MATRIX >::internal_get_block_vectors(), oomph::BlockPreconditioner< MATRIX >::internal_return_block_ordered_preconditioner_vector(), oomph::BlockPreconditioner< MATRIX >::internal_return_block_vector(), oomph::BlockPreconditioner< MATRIX >::internal_return_block_vectors(), oomph::PseudoElasticPreconditioner::lagrange_multiplier_preconditioner_solve(), oomph::PseudoElasticPreconditionerOld::lagrange_multiplier_preconditioner_solve(), oomph::CRDoubleMatrix::multiply(), oomph::DenseDoubleMatrix::multiply(), oomph::CCDoubleMatrix::multiply(), oomph::CRDoubleMatrix::multiply_transpose(), oomph::DenseDoubleMatrix::multiply_transpose(), oomph::CCDoubleMatrix::multiply_transpose(), oomph::Problem::newton_solve(), oomph::Problem::newton_solve_continuation(), operator+=(), operator-=(), operator==(), oomph::MatrixBasedDiagPreconditioner::preconditioner_solve(), oomph::LagrangeEnforcedFlowPreconditioner::preconditioner_solve(), oomph::GMRESBlockPreconditioner::preconditioner_solve(), oomph::DoubleMatrixBase::residual(), oomph::GMRES< MATRIX >::solve_helper(), oomph::AugmentedProblemGMRES::solve_helper(), oomph::DoubleVectorHelpers::split(), oomph::DoubleVectorHelpers::split_without_communication(), oomph::GMRESBlockPreconditioner::update(), oomph::GMRES< MATRIX >::update(), and oomph::AugmentedProblemGMRES::update().
|
inline |
access function to the underlying values (const version)
Definition at line 260 of file double_vector.h.
References Values_pt.
|
friend |
Ouput operator for DoubleVector.
Definition at line 949 of file double_vector.cc.
|
private |
indicates that the vector has been built and is usable
Definition at line 326 of file double_vector.h.
Referenced by build(), built(), clear(), initialise(), and set_external_values().
|
private |
Boolean flag to indicate whether the vector's data (values_pt) is owned by this vector.
Definition at line 323 of file double_vector.h.
Referenced by build(), clear(), redistribute(), and set_external_values().
|
private |
the local vector
Definition at line 319 of file double_vector.h.
Referenced by build(), clear(), dot(), initialise(), max(), norm(), operator*=(), operator+=(), operator-=(), operator==(), operator[](), output(), output_local_values(), output_local_values_with_offset(), redistribute(), set_external_values(), and values_pt().