Thyra Version of the Day
Loading...
Searching...
No Matches
ExampleTridiagSerialLinearOp.hpp
1// @HEADER
2// *****************************************************************************
3// Thyra: Interfaces and Support for Abstract Numerical Algorithms
4//
5// Copyright 2004 NTESS and the Thyra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef THYRA_EXAMPLE_TRIDIAG_SERIAL_LINEAR_OP_HPP
11#define THYRA_EXAMPLE_TRIDIAG_SERIAL_LINEAR_OP_HPP
12
13#include "Thyra_LinearOpDefaultBase.hpp"
14#include "Thyra_DefaultSpmdVectorSpace.hpp"
15#include "Thyra_DetachedVectorView.hpp"
16#include "Teuchos_Assert.hpp"
17
18
45template<class Scalar>
47{
48public:
49
52
55 const Thyra::Ordinal dim,
59 )
60 { this->initialize(dim, lower, diag, upper); }
61
84 const Thyra::Ordinal dim,
88 )
89 {
90 TEUCHOS_TEST_FOR_EXCEPT( dim < 2 );
91 space_ = Thyra::defaultSpmdVectorSpace<Scalar>(dim);
92 lower_ = lower;
93 diag_ = diag;
94 upper_ = upper;
95 }
96
97protected:
98
99 // Overridden from LinearOpBase
100
104
108
111 { return true; } // This class supports everything!
112
114 void applyImpl(
115 const Thyra::EOpTransp M_trans,
118 const Scalar alpha,
119 const Scalar beta
120 ) const;
121
122private:
123
125 Teuchos::Array<Scalar> lower_; // size = dim - 1
126 Teuchos::Array<Scalar> diag_; // size = dim
127 Teuchos::Array<Scalar> upper_; // size = dim - 1
128
129}; // end class ExampleTridiagSerialLinearOp
130
131
132template<class Scalar>
134 const Thyra::EOpTransp M_trans,
137 const Scalar alpha,
138 const Scalar beta
139 ) const
140{
141
143 typedef Thyra::Ordinal Ordinal;
144
145 const Ordinal dim = space_->dim();
146
147 // Loop over the input columns
148
149 const Ordinal m = X_in.domain()->dim();
150
151 for (Ordinal col_j = 0; col_j < m; ++col_j) {
152
153 // Get access the the elements of column j
155 Thyra::DetachedVectorView<Scalar> y_vec(Y_inout->col(col_j));
156 const Teuchos::ArrayRCP<const Scalar> x = x_vec.sv().values();
157 const Teuchos::ArrayRCP<Scalar> y = y_vec.sv().values();
158
159 // Perform y = beta*y (being careful to set y=0 if beta=0 in case y is
160 // uninitialized on input!)
161 if( beta == ST::zero() ) {
162 for( Ordinal k = 0; k < dim; ++k ) y[k] = ST::zero();
163 }
164 else if( beta != ST::one() ) {
165 for( Ordinal k = 0; k < dim; ++k ) y[k] *= beta;
166 }
167
168 // Perform y = alpha*op(M)*x
169 Ordinal k = 0;
170 if( M_trans == Thyra::NOTRANS ) {
171 y[k] += alpha * ( diag_[k]*x[k] + upper_[k]*x[k+1] ); // First row
172 for( k = 1; k < dim - 1; ++k ) // Middle rows
173 y[k] += alpha * ( lower_[k-1]*x[k-1] + diag_[k]*x[k] + upper_[k]*x[k+1] );
174 y[k] += alpha * ( lower_[k-1]*x[k-1] + diag_[k]*x[k] ); // Last row
175 }
176 else if( M_trans == Thyra::CONJ ) {
177 y[k] += alpha * ( ST::conjugate(diag_[k])*x[k] + ST::conjugate(upper_[k])*x[k+1] );
178 for( k = 1; k < dim - 1; ++k )
179 y[k] += alpha * ( ST::conjugate(lower_[k-1])*x[k-1]
180 + ST::conjugate(diag_[k])*x[k] + ST::conjugate(upper_[k])*x[k+1] );
181 y[k] += alpha * ( ST::conjugate(lower_[k-1])*x[k-1] + ST::conjugate(diag_[k])*x[k] );
182 }
183 else if( M_trans == Thyra::TRANS ) {
184 y[k] += alpha * ( diag_[k]*x[k] + lower_[k]*x[k+1] );
185 for( k = 1; k < dim - 1; ++k )
186 y[k] += alpha * ( upper_[k-1]*x[k-1] + diag_[k]*x[k] + lower_[k]*x[k+1] );
187 y[k] += alpha * ( upper_[k-1]*x[k-1] + diag_[k]*x[k] );
188 }
189 else if( M_trans == Thyra::CONJTRANS ) {
190 y[k] += alpha * ( ST::conjugate(diag_[k])*x[k] + ST::conjugate(lower_[k])*x[k+1] );
191 for( k = 1; k < dim - 1; ++k )
192 y[k] += alpha * ( ST::conjugate(upper_[k-1])*x[k-1]
193 + ST::conjugate(diag_[k])*x[k] + ST::conjugate(lower_[k])*x[k+1] );
194 y[k] += alpha * ( ST::conjugate(upper_[k-1])*x[k-1] + ST::conjugate(diag_[k])*x[k] );
195 }
196 else {
197 TEUCHOS_TEST_FOR_EXCEPT(true); // Throw exception if we get here!
198 }
199 }
200
201}
202
203
204#endif // THYRA_EXAMPLE_TRIDIAG_SERIAL_LINEAR_OP_HPP
Simple example subclass for serial tridiagonal matrices.
ExampleTridiagSerialLinearOp()
Construct to uninitialized.
bool opSupportedImpl(Thyra::EOpTransp M_trans) const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > domain() const
void initialize(const Thyra::Ordinal dim, const Teuchos::ArrayView< const Scalar > &lower, const Teuchos::ArrayView< const Scalar > &diag, const Teuchos::ArrayView< const Scalar > &upper)
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > range() const
ExampleTridiagSerialLinearOp(const Thyra::Ordinal dim, const Teuchos::ArrayView< const Scalar > &lower, const Teuchos::ArrayView< const Scalar > &diag, const Teuchos::ArrayView< const Scalar > &upper)
initialize().
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
const ArrayRCP< const Scalar > values() const
const ArrayRCP< Scalar > values() const
Create an explicit non-mutable (const) view of a VectorBase object.
const RTOpPack::ConstSubVectorView< Scalar > & sv() const
Returns the explicit view as an RTOpPack::ConstSubVectorView<Scalar> object.
Create an explicit mutable (non-const) view of a VectorBase object.
const RTOpPack::SubVectorView< Scalar > & sv() const
Returns the explicit view as an RTOpPack::ConstSubVectorView<Scalar> object.
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_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. `*.
@ TRANS
Use the transposed operator.
@ NOTRANS
Use the non-transposed operator.
@ CONJTRANS
Use the transposed operator with complex-conjugate clements (same as TRANS for real scalar types).
@ CONJ
Use the non-transposed operator with complex-conjugate elements (same as NOTRANS for real scalar type...