Thyra Version of the Day
Loading...
Searching...
No Matches
ExampleTridiagSerialLinearOp.hpp
1// @HEADER
2// ***********************************************************************
3//
4// Thyra: Interfaces and Support for Abstract Numerical Algorithms
5// Copyright (2004) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Roscoe A. Bartlett (bartlettra@ornl.gov)
38//
39// ***********************************************************************
40// @HEADER
41
42#ifndef THYRA_EXAMPLE_TRIDIAG_SERIAL_LINEAR_OP_HPP
43#define THYRA_EXAMPLE_TRIDIAG_SERIAL_LINEAR_OP_HPP
44
45#include "Thyra_LinearOpDefaultBase.hpp"
46#include "Thyra_DefaultSpmdVectorSpace.hpp"
47#include "Thyra_DetachedVectorView.hpp"
48#include "Teuchos_Assert.hpp"
49
50
78template<class Scalar>
80{
81public:
82
85
88 const Thyra::Ordinal dim,
92 )
93 { this->initialize(dim, lower, diag, upper); }
94
117 const Thyra::Ordinal dim,
121 )
122 {
123 TEUCHOS_TEST_FOR_EXCEPT( dim < 2 );
124 space_ = Thyra::defaultSpmdVectorSpace<Scalar>(dim);
125 lower_ = lower;
126 diag_ = diag;
127 upper_ = upper;
128 }
129
130protected:
131
132 // Overridden from LinearOpBase
133
137
141
144 { return true; } // This class supports everything!
145
147 void applyImpl(
148 const Thyra::EOpTransp M_trans,
151 const Scalar alpha,
152 const Scalar beta
153 ) const;
154
155private:
156
158 Teuchos::Array<Scalar> lower_; // size = dim - 1
159 Teuchos::Array<Scalar> diag_; // size = dim
160 Teuchos::Array<Scalar> upper_; // size = dim - 1
161
162}; // end class ExampleTridiagSerialLinearOp
163
164
165template<class Scalar>
167 const Thyra::EOpTransp M_trans,
170 const Scalar alpha,
171 const Scalar beta
172 ) const
173{
174
176 typedef Thyra::Ordinal Ordinal;
177
178 const Ordinal dim = space_->dim();
179
180 // Loop over the input columns
181
182 const Ordinal m = X_in.domain()->dim();
183
184 for (Ordinal col_j = 0; col_j < m; ++col_j) {
185
186 // Get access the the elements of column j
188 Thyra::DetachedVectorView<Scalar> y_vec(Y_inout->col(col_j));
189 const Teuchos::ArrayRCP<const Scalar> x = x_vec.sv().values();
190 const Teuchos::ArrayRCP<Scalar> y = y_vec.sv().values();
191
192 // Perform y = beta*y (being careful to set y=0 if beta=0 in case y is
193 // uninitialized on input!)
194 if( beta == ST::zero() ) {
195 for( Ordinal k = 0; k < dim; ++k ) y[k] = ST::zero();
196 }
197 else if( beta != ST::one() ) {
198 for( Ordinal k = 0; k < dim; ++k ) y[k] *= beta;
199 }
200
201 // Perform y = alpha*op(M)*x
202 Ordinal k = 0;
203 if( M_trans == Thyra::NOTRANS ) {
204 y[k] += alpha * ( diag_[k]*x[k] + upper_[k]*x[k+1] ); // First row
205 for( k = 1; k < dim - 1; ++k ) // Middle rows
206 y[k] += alpha * ( lower_[k-1]*x[k-1] + diag_[k]*x[k] + upper_[k]*x[k+1] );
207 y[k] += alpha * ( lower_[k-1]*x[k-1] + diag_[k]*x[k] ); // Last row
208 }
209 else if( M_trans == Thyra::CONJ ) {
210 y[k] += alpha * ( ST::conjugate(diag_[k])*x[k] + ST::conjugate(upper_[k])*x[k+1] );
211 for( k = 1; k < dim - 1; ++k )
212 y[k] += alpha * ( ST::conjugate(lower_[k-1])*x[k-1]
213 + ST::conjugate(diag_[k])*x[k] + ST::conjugate(upper_[k])*x[k+1] );
214 y[k] += alpha * ( ST::conjugate(lower_[k-1])*x[k-1] + ST::conjugate(diag_[k])*x[k] );
215 }
216 else if( M_trans == Thyra::TRANS ) {
217 y[k] += alpha * ( diag_[k]*x[k] + lower_[k]*x[k+1] );
218 for( k = 1; k < dim - 1; ++k )
219 y[k] += alpha * ( upper_[k-1]*x[k-1] + diag_[k]*x[k] + lower_[k]*x[k+1] );
220 y[k] += alpha * ( upper_[k-1]*x[k-1] + diag_[k]*x[k] );
221 }
222 else if( M_trans == Thyra::CONJTRANS ) {
223 y[k] += alpha * ( ST::conjugate(diag_[k])*x[k] + ST::conjugate(lower_[k])*x[k+1] );
224 for( k = 1; k < dim - 1; ++k )
225 y[k] += alpha * ( ST::conjugate(upper_[k-1])*x[k-1]
226 + ST::conjugate(diag_[k])*x[k] + ST::conjugate(lower_[k])*x[k+1] );
227 y[k] += alpha * ( ST::conjugate(upper_[k-1])*x[k-1] + ST::conjugate(diag_[k])*x[k] );
228 }
229 else {
230 TEUCHOS_TEST_FOR_EXCEPT(true); // Throw exception if we get here!
231 }
232 }
233
234}
235
236
237#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...