46 #ifndef ANASAZI_THYRA_DEBUG_ADAPTER_HPP 
   47 #define ANASAZI_THYRA_DEBUG_ADAPTER_HPP 
   55 #include <Thyra_DetachedMultiVectorView.hpp> 
   56 #include <Thyra_MultiVectorBase.hpp> 
   57 #include <Thyra_MultiVectorStdOps.hpp> 
   59 #include "Teuchos_Assert.hpp" 
   81   template<
class ScalarType>
 
   99       _timerCreate(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::create")),
 
  100       _timerClone(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::clone")),
 
  101       _timerDestroy(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::destroy")),
 
  102       _timerMvTimesMatAddMv(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvtimesmataddmv")),
 
  103       _timerMvTransMv(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvtransmv")),
 
  104       _timerMvAddMv(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvaddmv")),
 
  105       _timerMvDot(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvdot")),
 
  106       _timerMvNorm(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvnorm")),
 
  107       _timerMvScale(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvscale")),
 
  108       _timerSetBlock(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::setblock")),
 
  109       _timerMvInit(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvinit")),
 
  110       _timerMvRandom(Teuchos::TimeMonitor::getNewTimer(
"ThyraMultiVec::mvrandom"))
 
  156       std::vector<Teuchos::RCP<Teuchos::Time> >  myTimers = 
getTimers();
 
  167       std::vector<Teuchos::RCP<Teuchos::Time> >  myTimers = 
getTimers();
 
  181       std::vector<Teuchos::RCP<Teuchos::Time> >  myTimers = 
getTimers();
 
  195       std::vector<Teuchos::RCP<Teuchos::Time> >  myTimers = 
getTimers();
 
  209       std::vector<Teuchos::RCP<Teuchos::Time> >  myTimers = 
getTimers();
 
  255 #ifdef HAVE_ANASAZI_EXPERIMENTAL
 
  268 #ifdef HAVE_ANASAZI_EXPERIMENTAL
 
  330     std::vector<Teuchos::RCP<Teuchos::Time> > 
getTimers()
 const {
 
  331      std::vector<Teuchos::RCP<Teuchos::Time> > timers;
 
  332      timers.push_back( _timerCreate );
 
  333      timers.push_back( _timerClone );
 
  334      timers.push_back( _timerDestroy );
 
  335      timers.push_back( _timerMvTimesMatAddMv );
 
  336      timers.push_back( _timerMvTransMv );
 
  337      timers.push_back( _timerMvAddMv );
 
  338      timers.push_back( _timerMvDot );
 
  339      timers.push_back( _timerMvNorm );
 
  340      timers.push_back( _timerMvScale );
 
  341      timers.push_back( _timerSetBlock );
 
  342      timers.push_back( _timerMvInit );
 
  343      timers.push_back( _timerMvRandom );
 
  351      _timerCreate = timers[0];
 
  352      _timerClone = timers[1];
 
  353      _timerDestroy = timers[2];
 
  354      _timerMvTimesMatAddMv = timers[3];
 
  355      _timerMvTransMv = timers[4];
 
  356      _timerMvAddMv = timers[5];
 
  357      _timerMvDot = timers[6];
 
  358      _timerMvNorm = timers[7];
 
  359      _timerMvScale = timers[8];
 
  360      _timerSetBlock = timers[9];
 
  361      _timerMvInit = timers[10];
 
  362      _timerMvRandom = timers[11];
 
  395   template<
class ScalarType>
 
static void MvNorm(const TMVB &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec)
Compute the 2-norm of each individual vector of mv. Upon return, normvec[i] holds the value of ...
 
int GetNumberVecs() const 
Obtain the vector length of *this. 
 
const MultiVec< ScalarType > * CloneView(const std::vector< int > &index) const 
Creates a new ThyraMultiVec that shares the selected contents of *this. 
 
MultiVec< ScalarType > * CloneViewNonConst(const std::vector< int > &index)
Creates a new ThyraMultiVec that shares the selected contents of *this. 
 
static void MvDot(const TMVB &mv, const TMVB &A, std::vector< ScalarType > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
 
ThyraOp(const Teuchos::RCP< const Thyra::LinearOpBase< ScalarType > > &Op)
Basic constructor. Accepts reference-counted pointer to an Thyra_Operator. 
 
void Apply(const MultiVec< ScalarType > &X, MultiVec< ScalarType > &Y) const 
This method takes the Anasazi::MultiVec X and applies the operator to it resulting in the Anasazi::Mu...
 
Teuchos::RCP< Thyra::MultiVectorBase< ScalarType > > getRCP()
Return the reference-counted pointer held by this object. 
 
ptrdiff_t GetGlobalLength() const 
Obtain the number of vectors in *this. 
 
Virtual base class which defines basic traits for the operator type. 
 
Template specialization of Anasazi::OperatorTraits class using the Thyra::LinearOpBase virtual base c...
 
void MvDot(const MultiVec< ScalarType > &A, std::vector< ScalarType > &b) const 
Compute a vector b where the components are the individual dot-products, i.e.  where A[i] is the i-th...
 
void copyTimers(const std::vector< Teuchos::RCP< Teuchos::Time > > &timers)
Copy a std::vector<> of timers into this object. 
 
static void MvTimesMatAddMv(const ScalarType alpha, const TMVB &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, TMVB &mv)
Update mv with . 
 
void MvInit(ScalarType alpha)
Replace each element of the vectors in *this with alpha. 
 
Template specialization of Anasazi::MultiVecTraits class using the Thyra::MultiVectorBase class...
 
Interface for multivectors used by Anasazi' linear solvers. 
 
void MvNorm(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec) const 
Compute the 2-norm of each individual vector of *this. Upon return, normvec[i] holds the 2-norm of th...
 
void MvScale(const std::vector< ScalarType > &alpha)
Scale each element of the i-th vector in *this with alpha[i]. 
 
static void MvPrint(const TMVB &mv, std::ostream &os)
Print the mv multi-vector to the os output stream. 
 
ConjType
Enumerated types used to specify conjugation arguments. 
 
void MvRandom()
Fill the vectors in *this with random numbers. 
 
void MvTransMv(ScalarType alpha, const MultiVec< ScalarType > &A, Teuchos::SerialDenseMatrix< int, ScalarType > &B) const 
Compute a dense matrix B through the matrix-matrix multiply . 
 
static void MvTransMv(const ScalarType alpha, const TMVB &A, const TMVB &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
Compute a dense matrix B through the matrix-matrix multiply . 
 
MultiVec< ScalarType > * CloneCopy(const std::vector< int > &index) const 
Creates a new ThyraMultiVec and copies the selected contents of *this into the new vector (deep copy)...
 
Basic adapter class for Anasazi::MultiVec that uses Thyra::MultiVectorBase<ScalarType>. 
 
Traits class which defines basic operations on multivectors. 
 
ThyraMultiVec(const Teuchos::RCP< Thyra::MultiVectorBase< ScalarType > > &mv, std::vector< Teuchos::RCP< Teuchos::Time > > &timers)
Basic ThyraMultiVec constructor (wraps Thyra::MultiVectorBase<> object). 
 
Specializations of the Anasazi multi-vector and operator traits classes using Thyra base classes Line...
 
void MvPrint(std::ostream &os) const 
Print *this ThyraMultiVec. 
 
static void MvAddMv(const ScalarType alpha, const TMVB &A, const ScalarType beta, const TMVB &B, TMVB &mv)
Replace mv with . 
 
Anasazi header file which uses auto-configuration information to include necessary C++ headers...
 
static void MvRandom(TMVB &mv)
Replace the vectors in mv with random vectors. 
 
void MvAddMv(ScalarType alpha, const MultiVec< ScalarType > &A, ScalarType beta, const MultiVec< ScalarType > &B)
Replace *this with . 
 
Basic adapter class for Anasazi::Operator that uses Thyra_Operator. 
 
std::vector< Teuchos::RCP< Teuchos::Time > > getTimers() const 
Return a std::vector<> of timers held by this object. 
 
void MvScale(ScalarType alpha)
Scale each element of the vectors in *this with alpha. 
 
static int GetNumberVecs(const TMVB &mv)
Obtain the number of vectors in mv. 
 
Templated virtual class for creating operators that can interface with the Anasazi::OperatorTraits cl...
 
MultiVec< ScalarType > * CloneCopy() const 
Creates a new ThyraMultiVec and copies contents of *this into the new vector (deep copy)...
 
static void SetBlock(const TMVB &A, const std::vector< int > &index, TMVB &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index...
 
static Teuchos::RCP< TMVB > CloneViewNonConst(TMVB &mv, const std::vector< int > &index)
Creates a new MultiVectorBase that shares the selected contents of mv (shallow copy). 
 
static void MvInit(TMVB &mv, ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha. 
 
Types and exceptions used within Anasazi solvers and interfaces. 
 
ThyraMultiVec(const ThyraMultiVec< ScalarType > &mv)
Copy constructor. 
 
static Teuchos::RCP< TMVB > CloneCopy(const TMVB &mv)
Creates a new MultiVectorBase and copies contents of mv into the new vector (deep copy)...
 
virtual ~ThyraMultiVec()
Destructor. 
 
MultiVec< ScalarType > * Clone(const int numvecs) const 
Creates a new empty ThyraMultiVec containing numvecs columns. 
 
void SetBlock(const MultiVec< ScalarType > &A, const std::vector< int > &index)
Copy the vectors in A to a set of vectors in *this. 
 
void MvTimesMatAddMv(ScalarType alpha, const MultiVec< ScalarType > &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, ScalarType beta)
Update *this with . 
 
static ptrdiff_t GetGlobalLength(const TMVB &mv)
Obtain the vector length of mv. 
 
static void Apply(const Thyra::LinearOpBase< ScalarType > &Op, const Thyra::MultiVectorBase< ScalarType > &x, Thyra::MultiVectorBase< ScalarType > &y)
This method takes the MultiVectorBase x and applies the LinearOpBase Op to it resulting in the MultiV...
 
Interface for multivectors used by Anasazi's linear solvers. 
 
Teuchos::RCP< const Thyra::MultiVectorBase< ScalarType > > getRCP() const 
Return the const reference-counted pointer held by this object. 
 
static Teuchos::RCP< const TMVB > CloneView(const TMVB &mv, const std::vector< int > &index)
Creates a new const MultiVectorBase that shares the selected contents of mv (shallow copy)...
 
ThyraMultiVec(const Teuchos::RCP< Thyra::MultiVectorBase< ScalarType > > &mv)
Basic ThyraMultiVec constructor (wraps Thyra::MultiVectorBase<> object). 
 
Anasazi's templated virtual class for constructing an operator that can interface with the OperatorTr...
 
static void MvScale(TMVB &mv, const ScalarType alpha)
Scale each element of the vectors in *this with alpha. 
 
static Teuchos::RCP< TMVB > Clone(const TMVB &mv, const int numvecs)
Creates a new empty MultiVectorBase containing numvecs columns.