10#ifndef THYRA_EXAMPLE_TRIDIAG_SPMD_LINEAR_OP_HPP 
   11#define THYRA_EXAMPLE_TRIDIAG_SPMD_LINEAR_OP_HPP 
   13#include "Thyra_LinearOpDefaultBase.hpp" 
   14#include "Thyra_DefaultSpmdVectorSpace.hpp" 
   15#include "Thyra_DetachedSpmdVectorView.hpp" 
   16#include "Teuchos_DefaultComm.hpp" 
   17#include "Teuchos_CommHelpers.hpp" 
  106template<
class Scalar>
 
  121    { this->
initialize(comm, localDim, lower, diag, upper); }
 
 
  187  void communicate( 
const bool first, 
const bool last,
 
 
  195template<
class Scalar>
 
  206  space_ = Thyra::defaultSpmdVectorSpace<Scalar>(comm, localDim, -1);
 
 
  213template<
class Scalar>
 
  225  using Teuchos::outArg;
 
  230  const Scalar zero = ST::zero();
 
  231  const Ordinal localDim = space_->localSubDim();
 
  232  const Ordinal procRank = comm_->getRank();
 
  233  const Ordinal numProcs = comm_->getSize();
 
  234  const Ordinal m = X_in.
domain()->dim();
 
  238  for (Ordinal col_j = 0; col_j < m; ++col_j) {
 
  247    const bool first = (procRank == 0), last = ( procRank == numProcs-1 );
 
  251    communicate( first, last, x(), outArg(x_km1), outArg(x_kp1) );
 
  258      y[k] = alpha * ( (first?zero:lower_[lk]*x_km1) + diag_[k]*x[k]
 
  259        + upper_[k]*x[k+1] );
 
  262      for( k = 1; k < localDim - 1; ++lk, ++k )
 
  263        y[k] = alpha * ( lower_[lk]*x[k-1] + diag_[k]*x[k] + upper_[k]*x[k+1] );
 
  265      y[k] = alpha * ( lower_[lk]*x[k-1] + diag_[k]*x[k]
 
  266        + (last?zero:upper_[k]*x_kp1) );
 
  270      y[k] = alpha * ( (first?zero:lower_[lk]*x_km1) + diag_[k]*x[k]
 
  271        + upper_[k]*x[k+1] ) + beta*y[k];
 
  274      for( k = 1; k < localDim - 1; ++lk, ++k )
 
  275        y[k] = alpha * ( lower_[lk]*x[k-1] + diag_[k]*x[k] + upper_[k]*x[k+1] )
 
  278      y[k] = alpha * ( lower_[lk]*x[k-1] + diag_[k]*x[k]
 
  279        + (last?zero:upper_[k]*x_kp1) ) + beta*y[k];
 
 
  287template<
class Scalar>
 
  289  const bool first, 
const bool last,
 
  296  const Ordinal localDim = space_->localSubDim();
 
  297  const Ordinal procRank = comm_->getRank();
 
  298  const Ordinal numProcs = comm_->getSize();
 
  299  const bool isEven = (procRank % 2 == 0);
 
  300  const bool isOdd = !isEven;
 
  303  if(isEven && procRank+1 < numProcs) send(*comm_, x[localDim-1], procRank+1);
 
  304  if(isOdd && procRank-1 >= 0) receive(*comm_, procRank-1, &*x_km1);
 
  307  if( isOdd && procRank-1 >= 0 )         send(*comm_, x[0], procRank-1);
 
  308  if( isEven && procRank+1 < numProcs ) receive(*comm_, procRank+1, &*x_kp1);
 
  311  if (isOdd && procRank+1 < numProcs) send(*comm_, x[localDim-1], procRank+1);
 
  312  if (isEven && procRank-1 >= 0) receive(*comm_, procRank-1, &*x_km1);
 
  315  if (isEven && procRank-1 >= 0) send(*comm_, x[0], procRank-1);
 
  316  if (isOdd && procRank+1 < numProcs) receive(*comm_, procRank+1, &*x_kp1);
 
Simple example subclass for Spmd tridiagonal matrices.
 
bool opSupportedImpl(Thyra::EOpTransp M_trans) const
 
void initialize(const Teuchos::RCP< const Teuchos::Comm< Thyra::Ordinal > > &comm, const Thyra::Ordinal localDim, const Teuchos::ArrayView< const Scalar > &lower, const Teuchos::ArrayView< const Scalar > &diag, const Teuchos::ArrayView< const Scalar > &upper)
 
void applyImpl(const Thyra::EOpTransp M_trans, const Thyra::MultiVectorBase< Scalar > &X_in, const Teuchos::Ptr< Thyra::MultiVectorBase< Scalar > > &Y_inout, const Scalar alpha, const Scalar beta) const
 
ExampleTridiagSpmdLinearOp()
 
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > domain() const
 
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > range() const
 
ExampleTridiagSpmdLinearOp(const Teuchos::RCP< const Teuchos::Comm< Thyra::Ordinal > > &comm, const Thyra::Ordinal localDim, const Teuchos::ArrayView< const Scalar > &lower, const Teuchos::ArrayView< const Scalar > &diag, const Teuchos::ArrayView< const Scalar > &upper)
 
const ArrayRCP< const Scalar > values() const
 
const ArrayRCP< Scalar > values() const
 
Create an explicit detached non-mutable (const) view of all of the local elements on this process of ...
 
const RTOpPack::ConstSubVectorView< Scalar > & sv() const
 
Create an explicit detached mutable (non-const) view of all of the local elements on this process of ...
 
const RTOpPack::SubVectorView< Scalar > & sv() const
 
virtual RCP< const VectorSpaceBase< Scalar > > domain() const =0
Return a smart pointer for the domain space for this operator.
 
Node subclass that provides a good default implementation for the describe() function.
 
Interface for a collection of column vectors called a multi-vector.
 
RCP< const VectorBase< Scalar > > col(Ordinal j) const
Calls colImpl().
 
#define TEUCHOS_ASSERT(assertion_test)
 
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
 
EOpTransp
Enumeration for determining how a linear operator is applied. `*.
 
Teuchos::Ordinal Ordinal
Type for the dimension of a vector space. `*.
 
@ NOTRANS
Use the non-transposed operator.