Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_BlockedVector_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Xpetra: A linear algebra interface package
4 //
5 // Copyright 2012 NTESS and the Xpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef XPETRA_BLOCKEDVECTOR_DEF_HPP
11 #define XPETRA_BLOCKEDVECTOR_DEF_HPP
12 
14 
15 #include "Xpetra_BlockedMultiVector.hpp"
16 #include "Xpetra_Exceptions.hpp"
17 
18 namespace Xpetra {
19 
20 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
22  BlockedVector(const Teuchos::RCP<const Xpetra::BlockedMap<LocalOrdinal, GlobalOrdinal, Node>>& map, bool zeroOut)
23  : Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(map, 1, zeroOut) {}
24 
25 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
26 BlockedVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
27  BlockedVector(Teuchos::RCP<const Xpetra::BlockedMap<LocalOrdinal, GlobalOrdinal, Node>> bmap,
29  : Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(bmap, v) {}
30 
31 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35  : Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(mapExtractor, v) {}
36 
37 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40 
41 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
45  assign(rhs); // dispatch to protected virtual method
46  return *this;
47 }
48 
49 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
51  replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar& value) {
52  BlockedMultiVector::replaceGlobalValue(globalRow, vectorIndex, value);
53 }
54 
55 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57  sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar& value) {
58  BlockedMultiVector::sumIntoGlobalValue(globalRow, vectorIndex, value);
59 }
60 
61 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63  replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar& value) {
64  BlockedMultiVector::replaceLocalValue(myRow, vectorIndex, value);
65 }
66 
67 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
69  sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar& value) {
70  BlockedMultiVector::sumIntoLocalValue(myRow, vectorIndex, value);
71 }
72 
73 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75  replaceGlobalValue(GlobalOrdinal globalRow, const Scalar& value) {
76  BlockedMultiVector::replaceGlobalValue(globalRow, 0, value);
77 }
78 
79 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81  sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar& value) {
82  BlockedMultiVector::sumIntoGlobalValue(globalRow, 0, value);
83 }
84 
85 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87  replaceLocalValue(LocalOrdinal myRow, const Scalar& value) {
89 }
90 
91 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
93  sumIntoLocalValue(LocalOrdinal myRow, const Scalar& value) {
95 }
96 
97 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
99  putScalar(const Scalar& value) {
101 }
102 
103 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
104 Teuchos::RCP<const Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
106  getVector(size_t j) const {
108 }
109 
110 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
111 Teuchos::RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
115 }
116 
117 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
118 Teuchos::ArrayRCP<const Scalar>
120  getData(size_t j) const {
121  return BlockedMultiVector::getData(j);
122 }
123 
124 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
125 Teuchos::ArrayRCP<Scalar>
127  getDataNonConst(size_t j) {
129 }
130 
131 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
133  dot(const MultiVector& A, const Teuchos::ArrayView<Scalar>& dots) const {
134  BlockedMultiVector::dot(A, dots);
135  return;
136 }
137 
138 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139 Scalar
142  Teuchos::Array<Scalar> dots = Teuchos::Array<Scalar>(1);
143  BlockedMultiVector::dot(A, dots);
144  return dots[0];
145 }
146 
147 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
151  return;
152 }
153 
154 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
158  return;
159 }
160 
161 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
163  scale(const Scalar& alpha) {
165  return;
166 }
167 
168 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
170  scale(Teuchos::ArrayView<const Scalar> alpha) {
172  return;
173 }
174 
175 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
177  update(const Scalar& alpha,
179  const Scalar& beta) {
180  BlockedMultiVector::update(alpha, A, beta);
181  return;
182 }
183 
184 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
186  update(const Scalar& alpha,
188  const Scalar& beta,
190  const Scalar& gamma) {
191  BlockedMultiVector::update(alpha, A, beta, B, gamma);
192  return;
193 }
194 
195 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
198  norm1() const {
199  using Array = Teuchos::Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>;
200  Array norm = Array(1);
201  this->norm1(norm);
202  return norm[0];
203 }
204 
205 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
206 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
208  norm2() const {
209  Teuchos::Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> norm =
210  Teuchos::Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>(1);
211  this->norm2(norm);
212  return norm[0];
213 }
214 
215 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
216 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
218  normInf() const {
219  Teuchos::Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>
220  norm = Teuchos::Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>(1);
221  this->normInf(norm);
222  return norm[0];
223 }
224 
225 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
227  norm1(const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>& norms) const {
229 }
230 
231 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
233  norm2(const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>& norms) const {
235 }
236 
237 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
239  normInf(const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>& norms) const {
241 }
242 
243 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
245  meanValue(const Teuchos::ArrayView<Scalar>& /* means */) const {
246  throw Xpetra::Exceptions::RuntimeError("BlockedVector::meanValue: Not (yet) supported by BlockedVector.");
247 }
248 
249 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
250 Scalar
252  meanValue() const {
253  throw Xpetra::Exceptions::RuntimeError("BlockedVector::meanValue: Not (yet) supported by BlockedVector.");
254 }
255 
256 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258  multiply(Teuchos::ETransp /* transA */,
259  Teuchos::ETransp /* transB */,
260  const Scalar& /* alpha */,
263  const Scalar& /* beta */) {
264  throw Xpetra::Exceptions::RuntimeError("BlockedVector::multiply: Not (yet) supported by BlockedVector.");
265 }
266 
267 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
269  multiply(Teuchos::ETransp /* transA */,
270  Teuchos::ETransp /* transB */,
271  const Scalar& /* alpha */,
274  const Scalar& /* beta */) {
275  throw Xpetra::Exceptions::RuntimeError("BlockedVector::multiply: Not (yet) supported by BlockedVector.");
276 }
277 
278 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
280  elementWiseMultiply(Scalar /* scalarAB */,
283  Scalar /* scalarThis */) {
284  throw Xpetra::Exceptions::RuntimeError("BlockedVector::elementWiseMultiply: Not (yet) supported by BlockedVector.");
285 }
286 
287 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289  elementWiseMultiply(Scalar /* scalarAB */,
292  Scalar /* scalarThis */) {
293  XPETRA_TEST_FOR_EXCEPTION(B.getMap()->isSameAs(*(this->getMap())) == false,
295  "BlockedVector::elementWiseMultipy: B must have same blocked map than this.");
296  TEUCHOS_TEST_FOR_EXCEPTION(A.getMap()->getLocalNumElements() != B.getMap()->getLocalNumElements(),
298  "BlockedVector::elementWiseMultipy: A has "
299  << A.getMap()->getLocalNumElements() << " elements, B has " << B.getMap()->getLocalNumElements()
300  << ".");
301  TEUCHOS_TEST_FOR_EXCEPTION(A.getMap()->getGlobalNumElements() != B.getMap()->getGlobalNumElements(),
303  "BlockedVector::elementWiseMultipy: A has " << A.getMap()->getGlobalNumElements()
304  << " elements, B has "
305  << B.getMap()->getGlobalNumElements() << ".");
306 
307  RCP<const BlockedMap> bmap = this->getBlockedMap();
308  RCP<const Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> rcpA = Teuchos::rcpFromRef(A);
309  RCP<const Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> bmvec = Teuchos::rcpFromRef(B);
310  RCP<const BlockedVector> bbmvec = Teuchos::rcp_dynamic_cast<const BlockedVector>(bmvec);
311  TEUCHOS_TEST_FOR_EXCEPTION(bbmvec.is_null() == true,
313  "BlockedVector::elementWiseMultipy: B must be a BlockedVector.");
314 
315  // TODO implement me
316  /*RCP<Xpetra::MapExtractor<Scalar,LocalOrdinal,GlobalOrdinal,Node> > me = Teuchos::rcp(new
317  Xpetra::MapExtractor<Scalar,LocalOrdinal,GlobalOrdinal,Node>(bmap));
318 
319  for(size_t m = 0; m < bmap->getNumMaps(); m++) {
320  // TODO introduce BlockedVector objects and "skip" this expensive ExtractVector call
321  RCP<const Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > pd = me->ExtractVector(rcpA,m,bmap->getThyraMode());
322  XPETRA_TEST_FOR_EXCEPTION(pd->getMap()->isSameAs(*(this->getBlockedMap()->getMap(m,bmap->getThyraMode())))==false,
323  Xpetra::Exceptions::RuntimeError, "BlockedVector::elementWiseMultipy: sub map of B does not fit with sub map of this.");
324  this->getMultiVector(m,bmap->getThyraMode())->elementWiseMultiply(scalarAB,*pd,*(bbmvec->getMultiVector(m,bmap->getThyraMode())),scalarThis);
325  }*/
326 }
327 
328 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329 size_t
331  getNumVectors() const {
332  return 1;
333 }
334 
335 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336 size_t
340  "BlockedVector::getLocalLength: routine not implemented. It has no value as one must iterate on the partial vectors.");
341  TEUCHOS_UNREACHABLE_RETURN(0);
342 }
343 
344 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348  return this->getBlockedMap()->getFullMap()->getGlobalNumElements();
349 }
350 
351 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
355  "BlockedVector::isSameSize: routine not implemented. It has no value as one must iterate on the partial vectors.");
356  TEUCHOS_UNREACHABLE_RETURN(0);
357 }
358 
359 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
360 std::string
362  description() const {
363  return std::string("BlockedVector");
364 }
365 
366 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
368  describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel) const {
369  out << description() << std::endl;
370  for (size_t r = 0; r < this->getBlockedMap()->getNumMaps(); r++) {
371  getMultiVector(r)->describe(out, verbLevel);
372  }
373 }
374 
375 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
377  replaceMap(const RCP<const Map>& map) {
379 }
380 
381 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
384  const Import& /* importer */,
385  CombineMode /* CM */) {
386  throw Xpetra::Exceptions::RuntimeError("BlockedVector::doImport: Not supported by BlockedVector.");
387 }
388 
389 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
392  const Import& /* importer */,
393  CombineMode /* CM */) {
394  throw Xpetra::Exceptions::RuntimeError("BlockedVector::doExport: Not supported by BlockedVector.");
395 }
396 
397 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
400  const Export& /* exporter */,
401  CombineMode /* CM */) {
402  throw Xpetra::Exceptions::RuntimeError("BlockedVector::doImport: Not supported by BlockedVector.");
403 }
404 
405 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
408  const Export& /* exporter */,
409  CombineMode /* CM */) {
410  throw Xpetra::Exceptions::RuntimeError("BlockedVector::doExport: Not supported by BlockedVector.");
411 }
412 
413 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
415  setSeed(unsigned int seed) {
416  for (size_t r = 0; r < this->getBlockedMap()->getNumMaps(); ++r) {
417  getMultiVector(r)->setSeed(seed);
418  }
419 }
420 
421 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
423  randomize(bool bUseXpetraImplementation) {
424  for (size_t r = 0; r < this->getBlockedMap()->getNumMaps(); ++r) {
425  getMultiVector(r)->randomize(bUseXpetraImplementation);
426  }
427 }
428 
429 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
431  randomize(const Scalar& minVal, const Scalar& maxVal, bool bUseXpetraImplementation) {
432  for (size_t r = 0; r < this->getBlockedMap()->getNumMaps(); ++r) {
433  getMultiVector(r)->randomize(minVal, maxVal, bUseXpetraImplementation);
434  }
435 }
436 
437 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
440  {
442  }
443 }
444 
445 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
447  Xpetra_randomize(const Scalar& minVal, const Scalar& maxVal) {
448  {
450  }
451 }
452 
453 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
454 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
456  getMap() const {
457  XPETRA_MONITOR("BlockedVector::getMap");
458  return this->getBlockedMap();
459 }
460 
461 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
462 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
464  getMultiVector(size_t r) const {
466 }
467 
468 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
469 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
471  getMultiVector(size_t r, bool bThyraMode) const {
472  return BlockedMultiVector::getMultiVector(r, bThyraMode);
473 }
474 
475 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
479  bool bThyraMode) {
480  BlockedMultiVector::setMultiVector(r, v, bThyraMode);
481  return;
482 }
483 
484 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
485 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
487  Merge() const {
488  return BlockedMultiVector::Merge();
489 }
490 
491 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
495 }
496 
497 // template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
498 // virtual void BlockedVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
499 // assign (const XpetrA::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& rhs)
500 // {
501 // throw Xpetra::Exceptions::RuntimeError("BlockedVector::assign: Not supported by BlockedVector.");
502 // }
503 
504 } // namespace Xpetra
505 
506 #endif // XPETRA_BLOCKEDVECTOR_DEF_HPP
virtual Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
virtual void reciprocal(const MultiVector &)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
virtual Scalar meanValue() const
Compute mean (average) value of this Vector.
void setMultiVector(size_t r, Teuchos::RCP< const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> v, bool bThyraMode)
set partial Vector associated with block row r
virtual void setSeed(unsigned int seed)
Set seed for Random function.
virtual void multiply(Teuchos::ETransp, Teuchos::ETransp, const Scalar &, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Scalar &)
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
virtual void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
virtual void update(const Scalar &alpha, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
virtual Teuchos::RCP< const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
virtual void sumIntoGlobalValue(GlobalOrdinal, size_t, const Scalar &)
Add value to existing value, using global (row) index.
virtual Teuchos::RCP< const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
virtual void elementWiseMultiply(Scalar, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, Scalar)
Element-wise multiply of a Vector A with a MultiVector B.
virtual void scale(const Scalar &alpha)
Scale the current values of a vector, this = alpha*this.
virtual void replaceMap(const RCP< const Map > &map)
virtual size_t getNumVectors() const
Number of columns in the Vector.
virtual Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
virtual void replaceGlobalValue(GlobalOrdinal, size_t, const Scalar &)
Replace value, using global (row) index.
virtual std::string description() const
A simple one-line description of this object.
Exception throws to report errors in the internal logical of the program.
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
virtual void abs(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input vector in target: A = abs(this).
virtual void sumIntoLocalValue(LocalOrdinal, size_t, const Scalar &)
Add value to existing value, using local (row) index.
virtual void dot(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
virtual void update(const Scalar &alpha, const MultiVector &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
virtual void replaceLocalValue(LocalOrdinal, size_t, const Scalar &)
Replace value, using local (row) index.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
Teuchos::RCP< MultiVector > Merge() const
merge BlockedMultiVector blocks to a single MultiVector
virtual void reciprocal(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise reciprocal values of input vector in target, this(i,j) = 1/A(i,j).
virtual void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual bool isSameSize(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &) const
Local number of rows on the calling process.
virtual void assign(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
virtual void replaceMap(const RCP< const Map > &map)
Teuchos::RCP< const Map > getMap() const
Access function for the underlying Map this DistObject was constructed with.
virtual size_t getLocalLength() const
Local number of rows on the calling process.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
The Map describing the parallel distribution of this object.
virtual void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
virtual void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
virtual ~BlockedVector()
Destructor.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of vector.
void setMultiVector(size_t r, Teuchos::RCP< const MultiVector > v, bool bThyraMode)
set partial multivector associated with block row r
virtual void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
virtual void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
size_t global_size_t
Global size_t object.
virtual Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
virtual void abs(const MultiVector &)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
BlockedVector(const Teuchos::RCP< const BlockedMap > &map, bool zeroOut=true)
Constructor.
virtual void dot(const MultiVector &, const Teuchos::ArrayView< Scalar > &) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
Teuchos::RCP< MultiVector > getMultiVector(size_t r) const
return partial multivector associated with block row r
virtual void assign(const MultiVector &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Merge() const
merge BlockedVector blocks to a single Vector
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
virtual Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this vector.
#define XPETRA_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Compute 1-norm of vector.
virtual void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this vector.
virtual global_size_t getGlobalLength() const
Global number of rows in the Vector.
virtual void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
BlockedVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & operator=(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Assignment operator: Does a deep copy.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
virtual void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm in vector.
virtual void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Import &, CombineMode)
Import.
virtual void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Import &, CombineMode)
Export.
virtual void putScalar(const Scalar &value)
Set all values in the vector with the given value.
virtual void Xpetra_randomize()
Set vector values to random numbers. XPetra implementation.
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getMultiVector(size_t r) const
return partial Vector associated with block row r
virtual void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.