Xpetra Version of the Day
Loading...
Searching...
No Matches
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
18namespace Xpetra {
19
20template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
23 : Xpetra::BlockedMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(map, 1, zeroOut) {}
24
25template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
30
31template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36
37template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40
41template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
45 assign(rhs); // dispatch to protected virtual method
46 return *this;
47}
48
49template <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
55template <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
61template <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
67template <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
73template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75 replaceGlobalValue(GlobalOrdinal globalRow, const Scalar& value) {
76 BlockedMultiVector::replaceGlobalValue(globalRow, 0, value);
77}
78
79template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81 sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar& value) {
82 BlockedMultiVector::sumIntoGlobalValue(globalRow, 0, value);
83}
84
85template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87 replaceLocalValue(LocalOrdinal myRow, const Scalar& value) {
89}
90
91template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
93 sumIntoLocalValue(LocalOrdinal myRow, const Scalar& value) {
95}
96
97template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
102
103template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
109
110template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
116
117template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
123
124template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
130
131template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
137
138template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139Scalar
146
147template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
153
154template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
160
161template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
163 scale(const Scalar& alpha) {
165 return;
166}
167
168template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
174
175template <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
184template <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
195template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
204
205template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
214
215template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224
225template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
230
231template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236
237template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
242
243template <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
249template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
250Scalar
252 meanValue() const {
253 throw Xpetra::Exceptions::RuntimeError("BlockedVector::meanValue: Not (yet) supported by BlockedVector.");
254}
255
256template <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
267template <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
278template <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
287template <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();
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
328template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329size_t
334
335template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336size_t
338 getLocalLength() const {
340 "BlockedVector::getLocalLength: routine not implemented. It has no value as one must iterate on the partial vectors.");
342}
343
344template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
347 getGlobalLength() const {
348 return this->getBlockedMap()->getFullMap()->getGlobalNumElements();
349}
350
351template <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.");
357}
358
359template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
360std::string
362 description() const {
363 return std::string("BlockedVector");
364}
365
366template <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
375template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
380
381template <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
389template <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
397template <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
405template <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
413template <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
421template <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
429template <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
437template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
444
445template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
452
453template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
456 getMap() const {
457 XPETRA_MONITOR("BlockedVector::getMap");
458 return this->getBlockedMap();
459}
460
461template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
467
468template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
474
475template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
483
484template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
490
491template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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
#define XPETRA_MONITOR(funcName)
#define XPETRA_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool is_null() const
virtual void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
virtual void sumIntoLocalValue(LocalOrdinal, size_t, const Scalar &)
Add value to existing value, using local (row) index.
virtual Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
Teuchos::RCP< MultiVector > getMultiVector(size_t r) const
return partial multivector associated with block row r
void setMultiVector(size_t r, Teuchos::RCP< const MultiVector > v, bool bThyraMode)
set partial multivector associated with block row r
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 Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
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]).
virtual void reciprocal(const MultiVector &)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,...
virtual void sumIntoGlobalValue(GlobalOrdinal, size_t, const Scalar &)
Add value to existing value, using global (row) index.
virtual void replaceLocalValue(LocalOrdinal, size_t, const Scalar &)
Replace value, using local (row) index.
Teuchos::RCP< MultiVector > Merge() const
merge BlockedMultiVector blocks to a single MultiVector
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 void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual void assign(const MultiVector &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
virtual void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
virtual void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
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 abs(const MultiVector &)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
virtual void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
virtual void replaceMap(const RCP< const Map > &map)
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).
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Merge() const
merge BlockedVector blocks to a single Vector
virtual size_t getLocalLength() const
Local number of rows on the calling process.
virtual void replaceMap(const RCP< const Map > &map)
virtual std::string description() const
A simple one-line description of this object.
BlockedVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & operator=(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Assignment operator: Does a deep copy.
virtual void Xpetra_randomize()
Set vector values to random numbers. XPetra implementation.
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 Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this vector.
virtual void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
Teuchos::RCP< const Map > getMap() const
Access function for the underlying Map this DistObject was constructed with.
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 Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute 2-norm of vector.
virtual void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
virtual Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
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 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 putScalar(const Scalar &value)
Set all values in the vector with the given value.
virtual void setSeed(unsigned int seed)
Set seed for Random function.
virtual void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Import &, CombineMode)
Export.
BlockedVector(const Teuchos::RCP< const BlockedMap > &map, bool zeroOut=true)
Constructor.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Compute 1-norm of vector.
virtual void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
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.
virtual void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
virtual void scale(const Scalar &alpha)
Scale the current values of a vector, this = alpha*this.
virtual void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this vector.
virtual void elementWiseMultiply(Scalar, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, Scalar)
Multiply a Vector A elementwise by a MultiVector B.
virtual size_t getNumVectors() const
Number of columns in the Vector.
virtual void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &, const Import &, CombineMode)
Import.
virtual Scalar meanValue() const
Compute mean (average) value of this Vector.
virtual bool isSameSize(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &) const
Local number of rows on the calling process.
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 global_size_t getGlobalLength() const
Global number of rows in the Vector.
virtual void assign(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getMultiVector(size_t r) const
return partial Vector associated with block row r
virtual ~BlockedVector()
Destructor.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute Inf-norm in vector.
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 Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
The Map describing the parallel distribution 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.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.