00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 #ifndef XPETRA_TPETRAMULTIVECTOR_HPP
00047 #define XPETRA_TPETRAMULTIVECTOR_HPP
00048
00049
00050
00051 #include "Xpetra_TpetraConfigDefs.hpp"
00052
00053 #include "Xpetra_MultiVector.hpp"
00054
00055 #include "Xpetra_TpetraMap.hpp"
00056 #include "Xpetra_Utils.hpp"
00057 #include "Xpetra_TpetraImport.hpp"
00058 #include "Xpetra_TpetraExport.hpp"
00059
00060 #include "Tpetra_MultiVector.hpp"
00061 #include "Tpetra_Vector.hpp"
00062
00063 namespace Xpetra {
00064
00065
00066 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00067 const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> & toTpetra(const MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
00068
00069 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00070 Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> & toTpetra(MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
00071
00072 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00073
00074 template<class S, class LO, class GO, class N> class TpetraVector;
00075 #endif
00076
00077
00078
00079 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00080 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
00081
00082 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00083 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
00084
00085
00086 template <class Scalar = MultiVector<>::scalar_type,
00087 class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type,
00088 class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type,
00089 class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
00090 class TpetraMultiVector
00091 : public virtual MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >
00092 {
00093
00094
00095 typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraMultiVectorClass;
00096
00097 public:
00098
00100
00101
00103 TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
00104 : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), NumVectors, zeroOut))) { }
00105
00107 TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
00108 : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(Tpetra::createCopy(toTpetra(source))))) { }
00109
00111 TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors)
00112 : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), A, LDA, NumVectors))) { }
00113
00115 TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
00116 : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), ArrayOfPtrs, NumVectors))) { }
00117
00118
00120 virtual ~TpetraMultiVector() { }
00121
00123
00125
00126
00128 void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::replaceGlobalValue"); vec_->replaceGlobalValue(globalRow, vectorIndex, value); }
00129
00131 void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::sumIntoGlobalValue"); vec_->sumIntoGlobalValue(globalRow, vectorIndex, value); }
00132
00134 void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::replaceLocalValue"); vec_->replaceLocalValue(myRow, vectorIndex, value); }
00135
00137 void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::sumIntoLocalValue"); vec_->sumIntoLocalValue(myRow, vectorIndex, value); }
00138
00140 void putScalar(const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::putScalar"); vec_->putScalar(value); }
00141
00143 void reduce() { XPETRA_MONITOR("TpetraMultiVector::reduce"); vec_->reduce(); }
00144
00146
00148
00149
00151 Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const { XPETRA_MONITOR("TpetraMultiVector::getVector"); return toXpetra(vec_->getVector(j)); }
00152
00154 Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j) { XPETRA_MONITOR("TpetraMultiVector::getVectorNonConst"); return toXpetra(vec_->getVectorNonConst(j)); }
00155
00157 Teuchos::ArrayRCP< const Scalar > getData(size_t j) const { XPETRA_MONITOR("TpetraMultiVector::getData"); return vec_->getData(j); }
00158
00160 Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j) { XPETRA_MONITOR("TpetraMultiVector::getDataNonConst"); return vec_->getDataNonConst(j); }
00161
00163 void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const { XPETRA_MONITOR("TpetraMultiVector::get1dCopy"); vec_->get1dCopy(A, LDA); }
00164
00166 void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const { XPETRA_MONITOR("TpetraMultiVector::get2dCopy"); vec_->get2dCopy(ArrayOfPtrs); }
00167
00169 Teuchos::ArrayRCP< const Scalar > get1dView() const { XPETRA_MONITOR("TpetraMultiVector::get1dView"); return vec_->get1dView(); }
00170
00172 Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const { XPETRA_MONITOR("TpetraMultiVector::get2dView"); return vec_->get2dView(); }
00173
00175 Teuchos::ArrayRCP< Scalar > get1dViewNonConst() { XPETRA_MONITOR("TpetraMultiVector::get1dViewNonConst"); return vec_->get1dViewNonConst(); }
00176
00178 Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst() { XPETRA_MONITOR("TpetraMultiVector::get2dViewNonConst"); return vec_->get2dViewNonConst(); }
00179
00181
00183
00184
00186 void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const { XPETRA_MONITOR("TpetraMultiVector::dot"); vec_->dot(toTpetra(A), dots); }
00187
00189 void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::abs"); vec_->abs(toTpetra(A)); }
00190
00192 void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::reciprocal"); vec_->reciprocal(toTpetra(A)); }
00193
00195 void scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha); }
00196
00198 void scale(Teuchos::ArrayView< const Scalar > alpha) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha); }
00199
00201 void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha, toTpetra(A)); }
00202
00204 void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { XPETRA_MONITOR("TpetraMultiVector::update"); vec_->update(alpha, toTpetra(A), beta); }
00205
00207 void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma) { XPETRA_MONITOR("TpetraMultiVector::update"); vec_->update(alpha, toTpetra(A), beta, toTpetra(B), gamma); }
00208
00210 void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::norm1"); vec_->norm1(norms); }
00211
00213 void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::norm2"); vec_->norm2(norms); }
00214
00216 void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::normInf"); vec_->normInf(norms); }
00217
00219 void normWeighted(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &weights, const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::normWeighted"); vec_->normWeighted(toTpetra(weights), norms); }
00220
00222 void meanValue(const Teuchos::ArrayView< Scalar > &means) const { XPETRA_MONITOR("TpetraMultiVector::meanValue"); vec_->meanValue(means); }
00223
00225 void multiply(Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta) { XPETRA_MONITOR("TpetraMultiVector::multiply"); vec_->multiply(transA, transB, alpha, toTpetra(A), toTpetra(B), beta); }
00226
00228
00230
00231
00233 size_t getNumVectors() const { XPETRA_MONITOR("TpetraMultiVector::getNumVectors"); return vec_->getNumVectors(); }
00234
00236 size_t getLocalLength() const { XPETRA_MONITOR("TpetraMultiVector::getLocalLength"); return vec_->getLocalLength(); }
00237
00239 global_size_t getGlobalLength() const { XPETRA_MONITOR("TpetraMultiVector::getGlobalLength"); return vec_->getGlobalLength(); }
00240
00242
00244
00245
00247 std::string description() const { XPETRA_MONITOR("TpetraMultiVector::description"); return vec_->description(); }
00248
00250 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraMultiVector::describe"); vec_->describe(out, verbLevel); }
00251
00253
00255 void elementWiseMultiply(Scalar scalarAB, const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &B, Scalar scalarThis);
00256
00257
00259 void randomize(bool bUseXpetraImplementation = false) {
00260 XPETRA_MONITOR("TpetraMultiVector::randomize");
00261
00262 if(bUseXpetraImplementation)
00263 MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Xpetra_randomize();
00264 else
00265 vec_->randomize();
00266 }
00267
00268
00269
00270
00271 Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> > getMap() const { XPETRA_MONITOR("TpetraMultiVector::getMap"); return toXpetra(vec_->getMap()); }
00272
00273 void doImport(const DistObject< Scalar, LocalOrdinal,GlobalOrdinal,Node> &source, const Import<LocalOrdinal,GlobalOrdinal,Node> &importer, CombineMode CM) {
00274 XPETRA_MONITOR("TpetraMultiVector::doImport");
00275
00276 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments.");
00277 RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tSource.getTpetra_MultiVector();
00278 this->getTpetra_MultiVector()->doImport(*v, toTpetra(importer), toTpetra(CM));
00279 }
00280
00281 void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import<LocalOrdinal,GlobalOrdinal,Node>& importer, CombineMode CM) {
00282 XPETRA_MONITOR("TpetraMultiVector::doExport");
00283
00284 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments.");
00285 RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tDest.getTpetra_MultiVector();
00286 this->getTpetra_MultiVector()->doExport(*v, toTpetra(importer), toTpetra(CM));
00287
00288 }
00289
00290 void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) {
00291 XPETRA_MONITOR("TpetraMultiVector::doImport");
00292
00293 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments.");
00294 RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tSource.getTpetra_MultiVector();
00295 this->getTpetra_MultiVector()->doImport(*v, toTpetra(exporter), toTpetra(CM));
00296
00297 }
00298
00299 void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) {
00300 XPETRA_MONITOR("TpetraMultiVector::doExport");
00301
00302 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments.");
00303 RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tDest.getTpetra_MultiVector();
00304 this->getTpetra_MultiVector()->doExport(*v, toTpetra(exporter), toTpetra(CM));
00305
00306 }
00307
00308 void replaceMap(const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map) {
00309 this->getTpetra_MultiVector()->replaceMap(toTpetra(map));
00310 }
00311
00312 template<class Node2>
00313 RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > clone(const RCP<Node2> &node2) const {
00314 return RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >(new TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node2>(vec_->clone(node2)));
00315
00316 }
00317
00319
00321
00322
00324 TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) : vec_(vec) { }
00325
00327 RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_MultiVector() const { return vec_; }
00328
00330 void setSeed(unsigned int seed) { XPETRA_MONITOR("TpetraMultiVector::seedrandom"); Teuchos::ScalarTraits< Scalar >::seedrandom(seed); }
00331
00333
00334 protected:
00337 virtual void
00338 assign (const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& rhs)
00339 {
00340 typedef TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> this_type;
00341 const this_type* rhsPtr = dynamic_cast<const this_type*> (&rhs);
00342 TEUCHOS_TEST_FOR_EXCEPTION(
00343 rhsPtr == NULL, std::invalid_argument, "Xpetra::MultiVector::operator=:"
00344 " The left-hand side (LHS) of the assignment has a different type than "
00345 "the right-hand side (RHS). The LHS has type Xpetra::TpetraMultiVector"
00346 " (which means it wraps a Tpetra::MultiVector), but the RHS has some "
00347 "other type. This probably means that the RHS wraps an "
00348 "Epetra_MultiVector. Xpetra::MultiVector does not currently implement "
00349 "assignment from an Epetra object to a Tpetra object, though this could"
00350 " be added with sufficient interest.");
00351
00352 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TMV;
00353 RCP<const TMV> rhsImpl = rhsPtr->getTpetra_MultiVector ();
00354 RCP<TMV> lhsImpl = this->getTpetra_MultiVector ();
00355
00356 TEUCHOS_TEST_FOR_EXCEPTION(
00357 rhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
00358 "(in Xpetra::TpetraMultiVector::assign): *this (the right-hand side of "
00359 "the assignment) has a null RCP<Tpetra::MultiVector> inside. Please "
00360 "report this bug to the Xpetra developers.");
00361 TEUCHOS_TEST_FOR_EXCEPTION(
00362 lhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
00363 "(in Xpetra::TpetraMultiVector::assign): The left-hand side of the "
00364 "assignment has a null RCP<Tpetra::MultiVector> inside. Please report "
00365 "this bug to the Xpetra developers.");
00366
00367 Tpetra::deep_copy (*lhsImpl, *rhsImpl);
00368 }
00369
00370 private:
00372 RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > vec_;
00373
00374 };
00375
00376
00377 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00378 const Tpetra::MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(const MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
00379 typedef TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass;
00380 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, x, tX, "toTpetra");
00381 return *tX.getTpetra_MultiVector();
00382 }
00383
00384 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00385 Tpetra::MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
00386 typedef TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass;
00387 XPETRA_DYNAMIC_CAST( TpetraMultiVectorClass, x, tX, "toTpetra");
00388 return *tX.getTpetra_MultiVector();
00389 }
00390
00391
00392
00393
00394 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00395 RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
00396 if (!vec.is_null())
00397 return rcp(new TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node >(vec));
00398
00399 return Teuchos::null;
00400 }
00401
00402 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00403 RCP<const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
00404 if (!vec.is_null())
00405 return rcp(new TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node >(vec));
00406
00407 return Teuchos::null;
00408 }
00409
00410 }
00411
00412
00413
00414
00415 #include "Xpetra_TpetraVector.hpp"
00416
00417 namespace Xpetra {
00418 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00419 void TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::elementWiseMultiply(Scalar scalarAB, const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &B, Scalar scalarThis) {
00420 XPETRA_MONITOR("TpetraMultiVector::elementWiseMultiply");
00421
00422
00423
00424 typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> tpv;
00425 XPETRA_DYNAMIC_CAST(const tpv, A, tA, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
00426 XPETRA_DYNAMIC_CAST(const TpetraMultiVector, B, tB, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
00427 vec_->elementWiseMultiply(scalarAB, *tA.getTpetra_Vector(), *tB.getTpetra_MultiVector(), scalarThis);
00428 }
00429
00430 }
00431
00432 #define XPETRA_TPETRAMULTIVECTOR_SHORT
00433 #endif // XPETRA_TPETRAMULTIVECTOR_HPP