Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosXpetraAdapterMultiVector_MP_Vector.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef BELOS_XPETRA_ADAPTER_MULTIVECTOR_MP_VECTOR_HPP
11 #define BELOS_XPETRA_ADAPTER_MULTIVECTOR_MP_VECTOR_HPP
12 
13 #include "BelosXpetraAdapterMultiVector.hpp"
16 
17 #ifdef HAVE_XPETRA_TPETRA
18 
19 namespace Belos { // should be moved to Belos or Xpetra?
20 
21  using Teuchos::RCP;
22  using Teuchos::rcp;
23 
25  //
26  // Implementation of the Belos::MultiVecTraits for Xpetra::MultiVector.
27  //
29 
34  template<class Storage, class LO, class GO, class Node>
35  class MultiVecTraits<typename Storage::value_type,
36  Xpetra::MultiVector<Sacado::MP::Vector<Storage>,
37  LO,GO,Node> > {
38  public:
39  typedef typename Storage::ordinal_type s_ordinal;
40  typedef typename Storage::value_type BaseScalar;
42  typedef typename Tpetra::MultiVector<Scalar,LO,GO,Node>::dot_type dot_type;
43  typedef typename Tpetra::MultiVector<Scalar,LO,GO,Node>::mag_type mag_type;
44 
45  private:
46 
47  typedef Xpetra::TpetraMultiVector<Scalar,LO,GO,Node> TpetraMultiVector;
48  typedef MultiVecTraits<dot_type,Tpetra::MultiVector<Scalar,LO,GO,Node> > MultiVecTraitsTpetra;
49 
50  public:
51 
52 #ifdef HAVE_BELOS_XPETRA_TIMERS
53  static RCP<Teuchos::Time> mvTimesMatAddMvTimer_, mvTransMvTimer_;
54 #endif
55 
56  static RCP<Xpetra::MultiVector<Scalar,LO,GO,Node> > Clone( const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, const int numvecs )
57  {
58  if (mv.getMap()->lib() == Xpetra::UseTpetra)
59  return rcp(new TpetraMultiVector(MultiVecTraitsTpetra::Clone(toTpetra(mv), numvecs)));
60  }
61 
62  static RCP<Xpetra::MultiVector<Scalar,LO,GO,Node> > CloneCopy( const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
63  {
64  if (mv.getMap()->lib() == Xpetra::UseTpetra)
65  return rcp(new TpetraMultiVector(MultiVecTraitsTpetra::CloneCopy(toTpetra(mv))));
66  }
67 
68  static RCP<Xpetra::MultiVector<Scalar,LO,GO,Node> > CloneCopy( const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, const std::vector<int>& index )
69  {
70  if (mv.getMap()->lib() == Xpetra::UseTpetra)
71  return rcp(new TpetraMultiVector(MultiVecTraitsTpetra::CloneCopy(toTpetra(mv), index)));
72  }
73 
74  static RCP<Xpetra::MultiVector<Scalar,LO,GO,Node> >
75  CloneCopy (const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv,
76  const Teuchos::Range1D& index)
77  {
78  if (mv.getMap()->lib() == Xpetra::UseTpetra)
79  return rcp(new TpetraMultiVector(MultiVecTraitsTpetra::CloneCopy(toTpetra(mv), index)));
80  }
81 
82  static RCP<Xpetra::MultiVector<Scalar,LO,GO,Node> > CloneViewNonConst( Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, const std::vector<int>& index )
83  {
84  if (mv.getMap()->lib() == Xpetra::UseTpetra)
85  return rcp(new TpetraMultiVector(MultiVecTraitsTpetra::CloneViewNonConst(toTpetra(mv), index)));
86  }
87 
88  static RCP<Xpetra::MultiVector<Scalar,LO,GO,Node> >
89  CloneViewNonConst(Xpetra::MultiVector<Scalar,LO,GO,Node>& mv,
90  const Teuchos::Range1D& index)
91  {
92  if (mv.getMap()->lib() == Xpetra::UseTpetra)
93  return rcp(new TpetraMultiVector(MultiVecTraitsTpetra::CloneViewNonConst(toTpetra(mv), index)));
94  }
95 
96  static RCP<const Xpetra::MultiVector<Scalar,LO,GO,Node> > CloneView(const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, const std::vector<int>& index )
97  {
98  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
99  //TODO: double check if the const_cast is safe here.
100  RCP<const Tpetra::MultiVector<Scalar,LO,GO,Node> > r = MultiVecTraitsTpetra::CloneView(toTpetra(mv), index);
101  return rcp(new TpetraMultiVector(Teuchos::rcp_const_cast<Tpetra::MultiVector<Scalar,LO,GO,Node> >(r)));
102  }
103  }
104 
105  static RCP<const Xpetra::MultiVector<Scalar,LO,GO,Node> >
106  CloneView (const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv,
107  const Teuchos::Range1D& index)
108  {
109  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
110  //TODO: double check if the const_cast is safe here.
111  RCP<const Tpetra::MultiVector<Scalar,LO,GO,Node> > r = MultiVecTraitsTpetra::CloneView(toTpetra(mv), index);
112  return rcp(new TpetraMultiVector(Teuchos::rcp_const_cast<Tpetra::MultiVector<Scalar,LO,GO,Node> >(r)));
113  }
114  }
115 
116  static ptrdiff_t GetGlobalLength( const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
117  {
118  if (mv.getMap()->lib() == Xpetra::UseTpetra)
119  return MultiVecTraitsTpetra::GetGlobalLength(toTpetra(mv));
120  }
121 
122  static int GetNumberVecs( const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
123  {
124  if (mv.getMap()->lib() == Xpetra::UseTpetra)
125  return MultiVecTraitsTpetra::GetNumberVecs(toTpetra(mv));
126  }
127 
128  static bool HasConstantStride( const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
129  {
130  if (mv.getMap()->lib() == Xpetra::UseTpetra)
131  return MultiVecTraitsTpetra::HasConstantStride(toTpetra(mv));
132  }
133 
134  static void MvTimesMatAddMv( dot_type alpha, const Xpetra::MultiVector<Scalar,LO,GO,Node>& A,
136  dot_type beta, Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
137  {
138 #ifdef HAVE_BELOS_XPETRA_TIMERS
139  Teuchos::TimeMonitor lcltimer(*mvTimesMatAddMvTimer_);
140 #endif
141  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
142  MultiVecTraitsTpetra::MvTimesMatAddMv(alpha, toTpetra(A), B, beta, toTpetra(mv));
143  return;
144  }
145  }
146 
147  static void MvAddMv( Scalar alpha, const Xpetra::MultiVector<Scalar,LO,GO,Node>& A, Scalar beta, const Xpetra::MultiVector<Scalar,LO,GO,Node>& B, Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
148  {
149  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
150  MultiVecTraitsTpetra::MvAddMv(alpha, toTpetra(A), beta, toTpetra(B), toTpetra(mv));
151  return;
152  }
153  }
154 
155  static void MvScale ( Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, Scalar alpha )
156  {
157  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
158  MultiVecTraitsTpetra::MvScale(toTpetra(mv), alpha);
159  return;
160  }
161  }
162 
163  static void MvScale ( Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, const std::vector<BaseScalar>& alphas )
164  {
165  std::vector<Scalar> alphas_mp(alphas.size());
166  const size_t sz = alphas.size();
167  for (size_t i=0; i<sz; ++i)
168  alphas_mp[i] = alphas[i];
169  MvScale (mv, alphas_mp);
170  }
171 
172  static void MvScale ( Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, const std::vector<Scalar>& alphas )
173  {
174  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
175  MultiVecTraitsTpetra::MvScale(toTpetra(mv), alphas);
176  return;
177  }
178  }
179 
180  static void MvTransMv( dot_type alpha, const Xpetra::MultiVector<Scalar,LO,GO,Node>& A, const Xpetra::MultiVector<Scalar,LO,GO,Node>& B, Teuchos::SerialDenseMatrix<int,dot_type>& C)
181  {
182 #ifdef HAVE_BELOS_XPETRA_TIMERS
183  Teuchos::TimeMonitor lcltimer(*mvTransMvTimer_);
184 #endif
185 
186  if (A.getMap()->lib() == Xpetra::UseTpetra) {
187  MultiVecTraitsTpetra::MvTransMv(alpha, toTpetra(A), toTpetra(B), C);
188  return;
189  }
190  }
191 
192  static void MvDot( const Xpetra::MultiVector<Scalar,LO,GO,Node>& A, const Xpetra::MultiVector<Scalar,LO,GO,Node>& B, std::vector<dot_type> &dots)
193  {
194  if (A.getMap()->lib() == Xpetra::UseTpetra) {
195  MultiVecTraitsTpetra::MvDot(toTpetra(A), toTpetra(B), dots);
196  return;
197  }
198  }
199 
200  static void MvNorm(const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, std::vector<mag_type> &normvec, NormType type=TwoNorm)
201  {
202  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
203  MultiVecTraitsTpetra::MvNorm(toTpetra(mv), normvec, type);
204  return;
205  }
206  }
207 
208  static void SetBlock( const Xpetra::MultiVector<Scalar,LO,GO,Node>& A, const std::vector<int>& index, Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
209  {
210  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
211  MultiVecTraitsTpetra::SetBlock(toTpetra(A), index, toTpetra(mv));
212  return;
213  }
214  }
215 
216  static void
217  SetBlock (const Xpetra::MultiVector<Scalar,LO,GO,Node>& A,
218  const Teuchos::Range1D& index,
219  Xpetra::MultiVector<Scalar,LO,GO,Node>& mv)
220  {
221  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
222  MultiVecTraitsTpetra::SetBlock(toTpetra(A), index, toTpetra(mv));
223  return;
224  }
225  }
226 
227  static void
228  Assign (const Xpetra::MultiVector<Scalar,LO,GO,Node>& A,
229  Xpetra::MultiVector<Scalar,LO,GO,Node>& mv)
230  {
231  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
232  MultiVecTraitsTpetra::Assign(toTpetra(A), toTpetra(mv));
233  return;
234  }
235  }
236 
237  static void MvRandom( Xpetra::MultiVector<Scalar,LO,GO,Node>& mv )
238  {
239  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
240  MultiVecTraitsTpetra::MvRandom(toTpetra(mv));
241  return;
242  }
243  }
244 
245  static void MvInit( Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, Scalar alpha = Teuchos::ScalarTraits<Scalar>::zero() )
246  {
247  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
248  MultiVecTraitsTpetra::MvInit(toTpetra(mv), alpha);
249  return;
250  }
251  }
252 
253  static void MvPrint( const Xpetra::MultiVector<Scalar,LO,GO,Node>& mv, std::ostream& os )
254  {
255  if (mv.getMap()->lib() == Xpetra::UseTpetra) {
256  MultiVecTraitsTpetra::MvPrint(toTpetra(mv), os);
257  return;
258  }
259  }
260 
261  };
262 
263 } // end of Belos namespace
264 
265 #endif
266 
267 #endif
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)