|
Tpetra Matrix/Vector Services
Version of the Day
|
00001 // @HEADER 00002 // *********************************************************************** 00003 // 00004 // Tpetra: Templated Linear Algebra Services Package 00005 // Copyright (2008) Sandia Corporation 00006 // 00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00008 // the U.S. Government retains certain rights in this software. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00038 // 00039 // ************************************************************************ 00040 // @HEADER 00041 00042 #ifndef TPETRA_EXPERIMENTAL_BLOCKVIEW_HPP 00043 #define TPETRA_EXPERIMENTAL_BLOCKVIEW_HPP 00044 00047 00048 #include <Teuchos_ScalarTraits.hpp> 00049 #include <Teuchos_LAPACK.hpp> 00050 00051 00052 //#include "Teuchos_LAPACK_wrappers.hpp" 00053 //extern "C" {int DGETRF_F77(const int *, const int *, double *, const int*, int *, int*);} 00054 00055 namespace Tpetra { 00056 00068 namespace Experimental { 00069 00084 template<class Scalar, class LO> 00085 class LittleBlock { 00086 private: 00087 typedef Teuchos::ScalarTraits<Scalar> STS; 00088 00089 public: 00095 LittleBlock (Scalar* const A, 00096 const LO blockSize, 00097 const LO strideX, 00098 const LO strideY) : 00099 A_ (A), blockSize_ (blockSize), strideX_ (strideX), strideY_ (strideY) 00100 { 00101 } 00102 00104 LO getBlockSize () const { 00105 return blockSize_; 00106 } 00107 00109 Scalar* getRawPtr () const { 00110 return A_; 00111 } 00112 00114 Scalar& operator() (const LO i, const LO j) const { 00115 return A_[i * strideX_ + j * strideY_]; 00116 } 00117 00119 template<class LittleBlockType> 00120 void update (const Scalar& alpha, const LittleBlockType& X) const { 00121 for (LO j = 0; j < blockSize_; ++j) { 00122 for (LO i = 0; i < blockSize_; ++i) { 00123 (*this)(i,j) += alpha * X(i,j); 00124 } 00125 } 00126 } 00127 00129 template<class LittleBlockType> 00130 void assign (const LittleBlockType& X) const { 00131 for (LO j = 0; j < blockSize_; ++j) { 00132 for (LO i = 0; i < blockSize_; ++i) { 00133 (*this)(i,j) = X(i,j); 00134 } 00135 } 00136 } 00137 00139 void scale (const Scalar& alpha) const { 00140 for (LO j = 0; j < blockSize_; ++j) { 00141 for (LO i = 0; i < blockSize_; ++i) { 00142 (*this)(i,j) *= alpha; 00143 } 00144 } 00145 } 00146 00148 void fill (const Scalar& alpha) const { 00149 for (LO j = 0; j < blockSize_; ++j) { 00150 for (LO i = 0; i < blockSize_; ++i) { 00151 (*this)(i,j) = alpha; 00152 } 00153 } 00154 } 00155 00160 template<class LittleBlockType> 00161 void absmax (const LittleBlockType& X) const { 00162 for (LO j = 0; j < blockSize_; ++j) { 00163 for (LO i = 0; i < blockSize_; ++i) { 00164 Scalar& Y_ij = (*this)(i,j); 00165 const Scalar X_ij = X(i,j); 00166 Y_ij = std::max (STS::magnitude (Y_ij), STS::magnitude (X_ij)); 00167 } 00168 } 00169 } 00170 00171 00172 00173 void factorize (int* ipiv, int & info) 00174 { 00175 Teuchos::LAPACK<int, Scalar> lapack; 00176 lapack.GETRF(blockSize_, blockSize_, A_, blockSize_, ipiv, &info); 00177 } 00178 00179 template<class LittleVectorType> 00180 void solve (LittleVectorType & X, const int* ipiv) const 00181 { 00182 int info; 00183 Teuchos::LAPACK<int, Scalar> lapack; 00184 char trans = 'T'; 00185 lapack.GETRS(trans, blockSize_, 1, A_, blockSize_, ipiv, X.getRawPtr(), blockSize_, &info); 00186 } 00187 00188 private: 00189 Scalar* const A_; 00190 const LO blockSize_; 00191 const LO strideX_; 00192 const LO strideY_; 00193 std::vector<int> ipiv_; 00194 }; 00195 00196 00213 template<class Scalar, class LO> 00214 class LittleVector { 00215 private: 00216 typedef Teuchos::ScalarTraits<Scalar> STS; 00217 00218 public: 00223 LittleVector (Scalar* const A, const LO blockSize, const LO stride) : 00224 A_ (A), blockSize_ (blockSize), strideX_ (stride) 00225 {} 00226 00228 Scalar* getRawPtr () const { 00229 return A_; 00230 } 00231 00233 LO getBlockSize () const { 00234 return blockSize_; 00235 } 00236 00238 LO getStride () const { 00239 return strideX_; 00240 } 00241 00243 Scalar& operator() (const LO i) const { 00244 return A_[i * strideX_]; 00245 } 00246 00248 template<class LittleVectorType> 00249 void update (const Scalar& alpha, const LittleVectorType& X) const { 00250 for (LO i = 0; i < blockSize_; ++i) { 00251 (*this)(i) += alpha * X(i); 00252 } 00253 } 00254 00256 template<class LittleVectorType> 00257 void assign (const LittleVectorType& X) const { 00258 for (LO i = 0; i < blockSize_; ++i) { 00259 (*this)(i) = X(i); 00260 } 00261 } 00262 00264 void scale (const Scalar& alpha) const { 00265 for (LO i = 0; i < blockSize_; ++i) { 00266 (*this)(i) *= alpha; 00267 } 00268 } 00269 00271 void fill (const Scalar& alpha) const { 00272 for (LO i = 0; i < blockSize_; ++i) { 00273 (*this)(i) = alpha; 00274 } 00275 } 00276 00281 template<class LittleVectorType> 00282 void absmax (const LittleVectorType& X) const { 00283 for (LO i = 0; i < blockSize_; ++i) { 00284 Scalar& Y_i = (*this)(i); 00285 Y_i = std::max (STS::magnitude (Y_i), STS::magnitude (X (i))); 00286 } 00287 } 00288 00290 template<class LittleVectorType> 00291 bool equal (const LittleVectorType& X) const { 00292 if (getBlockSize () != X.getBlockSize ()) { 00293 return false; 00294 } 00295 for (LO i = 0; i < blockSize_; ++i) { 00296 if ((*this)(i) != X(i)) { 00297 return false; 00298 } 00299 } 00300 return true; 00301 } 00302 00304 template<class LittleBlockType, class LittleVectorType> 00305 void 00306 matvecUpdate (const Scalar& alpha, 00307 const LittleBlockType& A, 00308 const LittleVectorType& X) const 00309 { 00310 // FIXME (mfh 07 May 2014) This is suitable for column major, not 00311 // for row major. Of course, we'll have to change other loops 00312 // above as well to make row major faster. 00313 for (LO j = 0; j < blockSize_; ++j) { 00314 for (LO i = 0; i < blockSize_; ++i) { 00315 (*this)(i) += alpha * A(i,j) * X(j); 00316 } 00317 } 00318 } 00319 00320 private: 00321 Scalar* const A_; 00322 const LO blockSize_; 00323 const LO strideX_; 00324 }; 00325 00326 } // namespace Experimental 00327 } // namespace Tpetra 00328 00329 #endif // TPETRA_EXPERIMENTAL_BLOCKVIEW_HPP
1.7.6.1