/usr/include/trilinos/Xpetra_TpetraMultiVector.hpp is in libtrilinos-xpetra-dev 12.10.1-3.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 | // @HEADER
//
// ***********************************************************************
//
// Xpetra: A linear algebra interface package
// Copyright 2012 Sandia Corporation
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact
// Jonathan Hu (jhu@sandia.gov)
// Andrey Prokopenko (aprokop@sandia.gov)
// Ray Tuminaro (rstumin@sandia.gov)
//
// ***********************************************************************
//
// @HEADER
#ifndef XPETRA_TPETRAMULTIVECTOR_HPP
#define XPETRA_TPETRAMULTIVECTOR_HPP
/* this file is automatically generated - do not edit (see script/tpetra.py) */
#include "Xpetra_TpetraConfigDefs.hpp"
#include "Xpetra_MultiVector.hpp"
#include "Xpetra_TpetraMap.hpp" //TMP
#include "Xpetra_Utils.hpp"
#include "Xpetra_TpetraImport.hpp"
#include "Xpetra_TpetraExport.hpp"
#include "Tpetra_MultiVector.hpp"
#include "Tpetra_Vector.hpp"
namespace Xpetra {
// TODO: move that elsewhere
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> & toTpetra(const MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> & toTpetra(MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
#ifndef DOXYGEN_SHOULD_SKIP_THIS
// forward declaration of TpetraVector, needed to prevent circular inclusions
template<class S, class LO, class GO, class N> class TpetraVector;
#endif
// Because we aren't including the header...
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
template <class Scalar = MultiVector<>::scalar_type,
class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type,
class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type,
class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
class TpetraMultiVector
: public virtual MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >
{
// The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraMultiVectorClass;
public:
//! @name Constructors and destructor
//@{
//! Basic constuctor.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
: vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), NumVectors, zeroOut))) {
// TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,numVecs,zeroOut): numVecs = " << NumVectors << " < 1.");
}
//! Copy constructor (performs a deep copy).
TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
: vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(Tpetra::createCopy(toTpetra(source))))) { }
//! Create multivector by copying two-dimensional array of local data.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors)
: vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), A, LDA, NumVectors))) {
// TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,A,LDA,numVecs): numVecs = " << NumVectors << " < 1.");
}
//! Create multivector by copying array of views of local data.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
: vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), ArrayOfPtrs, NumVectors))) {
// TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,ArrayOfPtrs,numVecs): numVecs = " << NumVectors << " < 1.");
}
//! Destructor (virtual for memory safety of derived classes).
virtual ~TpetraMultiVector() { }
//@}
//! @name Post-construction modification routines
//@{
//! Replace value, using global (row) index.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::replaceGlobalValue"); vec_->replaceGlobalValue(globalRow, vectorIndex, value); }
//! Add value to existing value, using global (row) index.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::sumIntoGlobalValue"); vec_->sumIntoGlobalValue(globalRow, vectorIndex, value); }
//! Replace value, using local (row) index.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::replaceLocalValue"); vec_->replaceLocalValue(myRow, vectorIndex, value); }
//! Add value to existing value, using local (row) index.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::sumIntoLocalValue"); vec_->sumIntoLocalValue(myRow, vectorIndex, value); }
//! Set all values in the multivector with the given value.
void putScalar(const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::putScalar"); vec_->putScalar(value); }
//! Sum values of a locally replicated multivector across all processes.
void reduce() { XPETRA_MONITOR("TpetraMultiVector::reduce"); vec_->reduce(); }
//@}
//! @name Data Copy and View get methods
//@{
//! Return a Vector which is a const view of column j.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const { XPETRA_MONITOR("TpetraMultiVector::getVector"); return toXpetra(vec_->getVector(j)); }
//! Return a Vector which is a nonconst view of column j.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j) { XPETRA_MONITOR("TpetraMultiVector::getVectorNonConst"); return toXpetra(vec_->getVectorNonConst(j)); }
//! Const view of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const { XPETRA_MONITOR("TpetraMultiVector::getData"); return vec_->getData(j); }
//! View of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j) { XPETRA_MONITOR("TpetraMultiVector::getDataNonConst"); return vec_->getDataNonConst(j); }
//! Fill the given array with a copy of this multivector's local values.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const { XPETRA_MONITOR("TpetraMultiVector::get1dCopy"); vec_->get1dCopy(A, LDA); }
//! Fill the given array with a copy of this multivector's local values.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const { XPETRA_MONITOR("TpetraMultiVector::get2dCopy"); vec_->get2dCopy(ArrayOfPtrs); }
//! Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< const Scalar > get1dView() const { XPETRA_MONITOR("TpetraMultiVector::get1dView"); return vec_->get1dView(); }
//! Return const persisting pointers to values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const { XPETRA_MONITOR("TpetraMultiVector::get2dView"); return vec_->get2dView(); }
//! Nonconst persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst() { XPETRA_MONITOR("TpetraMultiVector::get1dViewNonConst"); return vec_->get1dViewNonConst(); }
//! Return non-const persisting pointers to values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst() { XPETRA_MONITOR("TpetraMultiVector::get2dViewNonConst"); return vec_->get2dViewNonConst(); }
//@}
//! @name Mathematical methods
//@{
//! Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const { XPETRA_MONITOR("TpetraMultiVector::dot"); vec_->dot(toTpetra(A), dots); }
//! Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::abs"); vec_->abs(toTpetra(A)); }
//! Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::reciprocal"); vec_->reciprocal(toTpetra(A)); }
//! Scale the current values of a multi-vector, this = alpha*this.
void scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha); }
//! Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void scale(Teuchos::ArrayView< const Scalar > alpha) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha); }
//! Replace multi-vector values with scaled values of A, this = alpha*A.
void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha, toTpetra(A)); }
//! Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
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); }
//! Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
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); }
//! Compute 1-norm of each vector in multi-vector.
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::norm1"); vec_->norm1(norms); }
//!
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::norm2"); vec_->norm2(norms); }
//! Compute Inf-norm of each vector in multi-vector.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::normInf"); vec_->normInf(norms); }
//! Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined for non-floating point scalar types (e.g., int).
void meanValue(const Teuchos::ArrayView< Scalar > &means) const { XPETRA_MONITOR("TpetraMultiVector::meanValue"); vec_->meanValue(means); }
//! Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
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); }
//@}
//! @name Attribute access functions
//@{
//! Number of columns in the multivector.
size_t getNumVectors() const { XPETRA_MONITOR("TpetraMultiVector::getNumVectors"); return vec_->getNumVectors(); }
//! Local number of rows on the calling process.
size_t getLocalLength() const { XPETRA_MONITOR("TpetraMultiVector::getLocalLength"); return vec_->getLocalLength(); }
//! Global number of rows in the multivector.
global_size_t getGlobalLength() const { XPETRA_MONITOR("TpetraMultiVector::getGlobalLength"); return vec_->getGlobalLength(); }
//@}
//! @name Overridden from Teuchos::Describable
//@{
//! A simple one-line description of this object.
std::string description() const { XPETRA_MONITOR("TpetraMultiVector::description"); return vec_->description(); }
//! Print the object with the given verbosity level to a FancyOStream.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraMultiVector::describe"); vec_->describe(out, verbLevel); }
//@}
//! Element-wise multiply of a Vector A with a TpetraMultiVector B.
void elementWiseMultiply(Scalar scalarAB, const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &B, Scalar scalarThis); // definition at the end of this file
//TODO: void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis){ vec_->elementWiseMultiply(scalarAB, toTpetra(A), toTpetra(B), scalarThis); }
//! Set multi-vector values to random numbers.
void randomize(bool bUseXpetraImplementation = false) {
XPETRA_MONITOR("TpetraMultiVector::randomize");
if(bUseXpetraImplementation)
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Xpetra_randomize();
else
vec_->randomize();
}
//{@
// Implements DistObject interface
Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> > getMap() const { XPETRA_MONITOR("TpetraMultiVector::getMap"); return toXpetra(vec_->getMap()); }
void doImport(const DistObject< Scalar, LocalOrdinal,GlobalOrdinal,Node> &source, const Import<LocalOrdinal,GlobalOrdinal,Node> &importer, CombineMode CM) {
XPETRA_MONITOR("TpetraMultiVector::doImport");
XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tSource.getTpetra_MultiVector();
this->getTpetra_MultiVector()->doImport(*v, toTpetra(importer), toTpetra(CM));
}
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import<LocalOrdinal,GlobalOrdinal,Node>& importer, CombineMode CM) {
XPETRA_MONITOR("TpetraMultiVector::doExport");
XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tDest.getTpetra_MultiVector();
this->getTpetra_MultiVector()->doExport(*v, toTpetra(importer), toTpetra(CM));
}
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) {
XPETRA_MONITOR("TpetraMultiVector::doImport");
XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tSource.getTpetra_MultiVector();
this->getTpetra_MultiVector()->doImport(*v, toTpetra(exporter), toTpetra(CM));
}
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) {
XPETRA_MONITOR("TpetraMultiVector::doExport");
XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal,Node> > v = tDest.getTpetra_MultiVector();
this->getTpetra_MultiVector()->doExport(*v, toTpetra(exporter), toTpetra(CM));
}
void replaceMap(const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map) {
XPETRA_MONITOR("TpetraMultiVector::replaceMap");
this->getTpetra_MultiVector()->replaceMap(toTpetra(map));
}
template<class Node2>
RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > clone(const RCP<Node2> &node2) const {
XPETRA_MONITOR("TpetraMultiVector::clone");
return RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >(new TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node2>(vec_->clone(node2)));
//toXpetra(vec_->clone(node2));
}
//@}
//! @name Xpetra specific
//@{
//! TpetraMultiVector constructor to wrap a Tpetra::MultiVector object
TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) : vec_(vec) { } //TODO removed const
//! Get the underlying Tpetra multivector
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_MultiVector() const { return vec_; }
//! Set seed for Random function.
void setSeed(unsigned int seed) { XPETRA_MONITOR("TpetraMultiVector::seedrandom"); Teuchos::ScalarTraits< Scalar >::seedrandom(seed); }
#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dual_view_type dual_view_type;
//typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::host_execution_space host_execution_space;
//typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dev_execution_space dev_execution_space;
/// \brief Return an unmanaged non-const view of the local data on a specific device.
/// \tparam TargetDeviceType The Kokkos Device type whose data to return.
///
/// \warning DO NOT USE THIS FUNCTION! There is no reason why you are working directly
/// with the Xpetra::TpetraMultiVector object. To write a code which is independent
/// from the underlying linear algebra package you should always use the abstract class,
/// i.e. Xpetra::MultiVector!
///
/// \warning Be aware that the view on the multivector data is non-persisting, i.e.
/// only valid as long as the multivector does not run of scope!
template<class TargetDeviceType>
typename Kokkos::Impl::if_c<
Kokkos::Impl::is_same<
typename dual_view_type::t_dev_um::execution_space::memory_space,
typename TargetDeviceType::memory_space>::value,
typename dual_view_type::t_dev_um,
typename dual_view_type::t_host_um>::type
getLocalView () const {
return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
}
typename dual_view_type::t_host_um getHostLocalView () const {
return subview(vec_->template getLocalView<typename dual_view_type::host_mirror_space> (),
Kokkos::ALL(), Kokkos::ALL());
}
typename dual_view_type::t_dev_um getDeviceLocalView() const {
return subview(vec_->template getLocalView<typename dual_view_type::t_dev_um::execution_space> (),
Kokkos::ALL(), Kokkos::ALL());
}
#endif
//@}
protected:
/// \brief Implementation of the assignment operator (operator=);
/// does a deep copy.
virtual void
assign (const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& rhs)
{
typedef TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> this_type;
const this_type* rhsPtr = dynamic_cast<const this_type*> (&rhs);
TEUCHOS_TEST_FOR_EXCEPTION(
rhsPtr == NULL, std::invalid_argument, "Xpetra::MultiVector::operator=:"
" The left-hand side (LHS) of the assignment has a different type than "
"the right-hand side (RHS). The LHS has type Xpetra::TpetraMultiVector"
" (which means it wraps a Tpetra::MultiVector), but the RHS has some "
"other type. This probably means that the RHS wraps an "
"Epetra_MultiVector. Xpetra::MultiVector does not currently implement "
"assignment from an Epetra object to a Tpetra object, though this could"
" be added with sufficient interest.");
typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TMV;
RCP<const TMV> rhsImpl = rhsPtr->getTpetra_MultiVector ();
RCP<TMV> lhsImpl = this->getTpetra_MultiVector ();
TEUCHOS_TEST_FOR_EXCEPTION(
rhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
"(in Xpetra::TpetraMultiVector::assign): *this (the right-hand side of "
"the assignment) has a null RCP<Tpetra::MultiVector> inside. Please "
"report this bug to the Xpetra developers.");
TEUCHOS_TEST_FOR_EXCEPTION(
lhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
"(in Xpetra::TpetraMultiVector::assign): The left-hand side of the "
"assignment has a null RCP<Tpetra::MultiVector> inside. Please report "
"this bug to the Xpetra developers.");
Tpetra::deep_copy (*lhsImpl, *rhsImpl);
}
private:
//! The Tpetra::MultiVector which this class wraps.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > vec_;
}; // TpetraMultiVector class
// TODO: move that elsewhere
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
const Tpetra::MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(const MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
typedef TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass;
XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, x, tX, "toTpetra");
return *tX.getTpetra_MultiVector();
}
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(MultiVector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
typedef TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass;
XPETRA_DYNAMIC_CAST( TpetraMultiVectorClass, x, tX, "toTpetra");
return *tX.getTpetra_MultiVector();
}
//
// Things we actually need
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
if (!vec.is_null())
return rcp(new TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node >(vec));
return Teuchos::null;
}
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP<const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
if (!vec.is_null())
return rcp(new TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node >(vec));
return Teuchos::null;
}
#ifdef HAVE_XPETRA_EPETRA
#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
(!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
// specialization for TpetraMultiVector on EpetraNode and GO=int
template <class Scalar>
class TpetraMultiVector<Scalar,int,int,EpetraNode>
: public virtual MultiVector< Scalar, int, int, EpetraNode >
{
typedef int LocalOrdinal;
typedef int GlobalOrdinal;
typedef EpetraNode Node;
// The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraMultiVectorClass;
public:
//! @name Constructors and destructor
//@{
//! Default constructor
TpetraMultiVector () {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "int", typeid(EpetraNode).name() );
}
//! Basic constuctor.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "int", typeid(EpetraNode).name() );
}
//! Copy constructor (performs a deep copy).
TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "int", typeid(EpetraNode).name() );
}
//! Create multivector by copying two-dimensional array of local data.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "int", typeid(EpetraNode).name() );
}
//! Create multivector by copying array of views of local data.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "int", typeid(EpetraNode).name() );
}
//! Destructor (virtual for memory safety of derived classes).
virtual ~TpetraMultiVector() { }
//@}
//! @name Post-construction modification routines
//@{
//! Replace value, using global (row) index.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
//! Add value to existing value, using global (row) index.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
//! Replace value, using local (row) index.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
//! Add value to existing value, using local (row) index.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
//! Set all values in the multivector with the given value.
void putScalar(const Scalar &value) { }
//! Sum values of a locally replicated multivector across all processes.
void reduce() { }
//@}
//! @name Data Copy and View get methods
//@{
//! Return a Vector which is a const view of column j.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const { return Teuchos::null; }
//! Return a Vector which is a nonconst view of column j.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j) { return Teuchos::null; }
//! Const view of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const { return Teuchos::ArrayRCP< const Scalar >(); }
//! View of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j) { return Teuchos::ArrayRCP< Scalar >(); }
//! Fill the given array with a copy of this multivector's local values.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const { }
//! Fill the given array with a copy of this multivector's local values.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const { }
//! Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< const Scalar > get1dView() const { return Teuchos::ArrayRCP< const Scalar >(); }
//! Return const persisting pointers to values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const { return Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > >(); }
//! Nonconst persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst() { return Teuchos::ArrayRCP< Scalar >(); }
//! Return non-const persisting pointers to values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst() { return Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > >(); }
//@}
//! @name Mathematical methods
//@{
//! Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const { }
//! Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
//! Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
//! Scale the current values of a multi-vector, this = alpha*this.
void scale(const Scalar &alpha) { }
//! Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void scale(Teuchos::ArrayView< const Scalar > alpha) { }
//! Replace multi-vector values with scaled values of A, this = alpha*A.
void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
//! Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { }
//! Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
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) { }
//! Compute 1-norm of each vector in multi-vector.
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
//!
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
//! Compute Inf-norm of each vector in multi-vector.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
//! Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined for non-floating point scalar types (e.g., int).
void meanValue(const Teuchos::ArrayView< Scalar > &means) const { }
//! Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
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) { }
//@}
//! @name Attribute access functions
//@{
//! Number of columns in the multivector.
size_t getNumVectors() const { return 0; }
//! Local number of rows on the calling process.
size_t getLocalLength() const { return 0; }
//! Global number of rows in the multivector.
global_size_t getGlobalLength() const { return 0; }
//@}
//! @name Overridden from Teuchos::Describable
//@{
//! A simple one-line description of this object.
std::string description() const { return std::string(""); }
//! Print the object with the given verbosity level to a FancyOStream.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { }
//@}
//! Element-wise multiply of a Vector A with a TpetraMultiVector B.
void elementWiseMultiply(Scalar scalarAB, const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &B, Scalar scalarThis) {};
//! Set multi-vector values to random numbers.
void randomize(bool bUseXpetraImplementation = false) { }
//{@
// Implements DistObject interface
Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> > getMap() const { return Teuchos::null; }
void doImport(const DistObject< Scalar, LocalOrdinal,GlobalOrdinal,Node> &source, const Import<LocalOrdinal,GlobalOrdinal,Node> &importer, CombineMode CM) { }
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import<LocalOrdinal,GlobalOrdinal,Node>& importer, CombineMode CM) { }
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) { }
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) { }
void replaceMap(const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map) { }
template<class Node2>
RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > clone(const RCP<Node2> &node2) const { return Teuchos::null; }
//@}
//! @name Xpetra specific
//@{
//! TpetraMultiVector constructor to wrap a Tpetra::MultiVector object
TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "int", typeid(EpetraNode).name() );
}
//! Get the underlying Tpetra multivector
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_MultiVector() const { return Teuchos::null; }
//! Set seed for Random function.
void setSeed(unsigned int seed) { }
#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dual_view_type dual_view_type;
//typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::host_execution_space host_execution_space;
//typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dev_execution_space dev_execution_space;
/// \brief Return an unmanaged non-const view of the local data on a specific device.
/// \tparam TargetDeviceType The Kokkos Device type whose data to return.
///
/// \warning DO NOT USE THIS FUNCTION! There is no reason why you are working directly
/// with the Xpetra::TpetraMultiVector object. To write a code which is independent
/// from the underlying linear algebra package you should always use the abstract class,
/// i.e. Xpetra::MultiVector!
///
/// \warning Be aware that the view on the multivector data is non-persisting, i.e.
/// only valid as long as the multivector does not run of scope!
template<class TargetDeviceType>
typename Kokkos::Impl::if_c<
Kokkos::Impl::is_same<
typename dual_view_type::t_dev_um::execution_space::memory_space,
typename TargetDeviceType::memory_space>::value,
typename dual_view_type::t_dev_um,
typename dual_view_type::t_host_um>::type
getLocalView () const {
typename Kokkos::Impl::if_c<
Kokkos::Impl::is_same<
typename dual_view_type::t_dev_um::execution_space::memory_space,
typename TargetDeviceType::memory_space>::value,
typename dual_view_type::t_dev_um,
typename dual_view_type::t_host_um>::type dummy;
return dummy;
}
typename dual_view_type::t_host_um getHostLocalView () const {
//return subview(vec_->template getLocalView<typename dual_view_type::host_mirror_space> (),
// Kokkos::ALL(), Kokkos::ALL());
return typename dual_view_type::t_host_um();
}
typename dual_view_type::t_dev_um getDeviceLocalView() const {
//return subview(vec_->template getLocalView<typename dual_view_type::t_dev_um::execution_space> (),
// Kokkos::ALL(), Kokkos::ALL());
return typename dual_view_type::t_dev_um();
}
#endif
//@}
protected:
/// \brief Implementation of the assignment operator (operator=);
/// does a deep copy.
virtual void
assign (const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& rhs)
{ }
}; // TpetraMultiVector class (specialization GO=int, NO=EpetraNode)
#endif
#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
(!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
// specialization for TpetraMultiVector on EpetraNode and GO=long long
template <class Scalar>
class TpetraMultiVector<Scalar,int,long long,EpetraNode>
: public virtual MultiVector< Scalar, int, long long, EpetraNode >
{
typedef int LocalOrdinal;
typedef long long GlobalOrdinal;
typedef EpetraNode Node;
// The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraMultiVectorClass;
public:
//! @name Constructors and destructor
//@{
//! Basic constuctor.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "long long", typeid(EpetraNode).name() );
}
//! Copy constructor (performs a deep copy).
TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "long long", typeid(EpetraNode).name() );
}
//! Create multivector by copying two-dimensional array of local data.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "long long", typeid(EpetraNode).name() );
}
//! Create multivector by copying array of views of local data.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "long long", typeid(EpetraNode).name() );
}
//! Destructor (virtual for memory safety of derived classes).
virtual ~TpetraMultiVector() { }
//@}
//! @name Post-construction modification routines
//@{
//! Replace value, using global (row) index.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
//! Add value to existing value, using global (row) index.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
//! Replace value, using local (row) index.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
//! Add value to existing value, using local (row) index.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
//! Set all values in the multivector with the given value.
void putScalar(const Scalar &value) { }
//! Sum values of a locally replicated multivector across all processes.
void reduce() { }
//@}
//! @name Data Copy and View get methods
//@{
//! Return a Vector which is a const view of column j.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const { return Teuchos::null; }
//! Return a Vector which is a nonconst view of column j.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j) { return Teuchos::null; }
//! Const view of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const { return Teuchos::ArrayRCP< const Scalar >(); }
//! View of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j) { return Teuchos::ArrayRCP< Scalar >(); }
//! Fill the given array with a copy of this multivector's local values.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const { }
//! Fill the given array with a copy of this multivector's local values.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const { }
//! Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< const Scalar > get1dView() const { return Teuchos::ArrayRCP< const Scalar >(); }
//! Return const persisting pointers to values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const { return Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > >(); }
//! Nonconst persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst() { return Teuchos::ArrayRCP< Scalar >(); }
//! Return non-const persisting pointers to values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst() { return Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > >(); }
//@}
//! @name Mathematical methods
//@{
//! Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const { }
//! Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
//! Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
//! Scale the current values of a multi-vector, this = alpha*this.
void scale(const Scalar &alpha) { }
//! Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void scale(Teuchos::ArrayView< const Scalar > alpha) { }
//! Replace multi-vector values with scaled values of A, this = alpha*A.
void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
//! Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { }
//! Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
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) { }
//! Compute 1-norm of each vector in multi-vector.
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
//!
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
//! Compute Inf-norm of each vector in multi-vector.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
//! Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined for non-floating point scalar types (e.g., int).
void meanValue(const Teuchos::ArrayView< Scalar > &means) const { }
//! Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
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) { }
//@}
//! @name Attribute access functions
//@{
//! Number of columns in the multivector.
size_t getNumVectors() const { return 0; }
//! Local number of rows on the calling process.
size_t getLocalLength() const { return 0; }
//! Global number of rows in the multivector.
global_size_t getGlobalLength() const { return 0; }
//@}
//! @name Overridden from Teuchos::Describable
//@{
//! A simple one-line description of this object.
std::string description() const { return std::string(""); }
//! Print the object with the given verbosity level to a FancyOStream.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { }
//@}
//! Element-wise multiply of a Vector A with a TpetraMultiVector B.
void elementWiseMultiply(Scalar scalarAB, const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &B, Scalar scalarThis) {};
//! Set multi-vector values to random numbers.
void randomize(bool bUseXpetraImplementation = false) { }
//{@
// Implements DistObject interface
Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> > getMap() const { return Teuchos::null; }
void doImport(const DistObject< Scalar, LocalOrdinal,GlobalOrdinal,Node> &source, const Import<LocalOrdinal,GlobalOrdinal,Node> &importer, CombineMode CM) { }
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import<LocalOrdinal,GlobalOrdinal,Node>& importer, CombineMode CM) { }
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) { }
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export<LocalOrdinal,GlobalOrdinal,Node>& exporter, CombineMode CM) { }
void replaceMap(const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map) { }
template<class Node2>
RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > clone(const RCP<Node2> &node2) const { return Teuchos::null; }
//@}
//! @name Xpetra specific
//@{
//! TpetraMultiVector constructor to wrap a Tpetra::MultiVector object
TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) {
XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "long long", typeid(EpetraNode).name() );
}
//! Get the underlying Tpetra multivector
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_MultiVector() const { return Teuchos::null; }
//! Set seed for Random function.
void setSeed(unsigned int seed) { }
#ifdef HAVE_XPETRA_KOKKOS_REFACTOR
typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dual_view_type dual_view_type;
//typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::host_execution_space host_execution_space;
//typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dev_execution_space dev_execution_space;
/// \brief Return an unmanaged non-const view of the local data on a specific device.
/// \tparam TargetDeviceType The Kokkos Device type whose data to return.
///
/// \warning DO NOT USE THIS FUNCTION! There is no reason why you are working directly
/// with the Xpetra::TpetraMultiVector object. To write a code which is independent
/// from the underlying linear algebra package you should always use the abstract class,
/// i.e. Xpetra::MultiVector!
///
/// \warning Be aware that the view on the multivector data is non-persisting, i.e.
/// only valid as long as the multivector does not run of scope!
template<class TargetDeviceType>
typename Kokkos::Impl::if_c<
Kokkos::Impl::is_same<
typename dual_view_type::t_dev_um::execution_space::memory_space,
typename TargetDeviceType::memory_space>::value,
typename dual_view_type::t_dev_um,
typename dual_view_type::t_host_um>::type
getLocalView () const {
typename Kokkos::Impl::if_c<
Kokkos::Impl::is_same<
typename dual_view_type::t_dev_um::execution_space::memory_space,
typename TargetDeviceType::memory_space>::value,
typename dual_view_type::t_dev_um,
typename dual_view_type::t_host_um>::type dummy;
return dummy;
}
typename dual_view_type::t_host_um getHostLocalView () const {
//return subview(vec_->template getLocalView<typename dual_view_type::host_mirror_space> (),
// Kokkos::ALL(), Kokkos::ALL());
return typename dual_view_type::t_host_um();
}
typename dual_view_type::t_dev_um getDeviceLocalView() const {
//return subview(vec_->template getLocalView<typename dual_view_type::t_dev_um::execution_space> (),
// Kokkos::ALL(), Kokkos::ALL());
return typename dual_view_type::t_dev_um();
}
#endif
//@}
protected:
/// \brief Implementation of the assignment operator (operator=);
/// does a deep copy.
virtual void
assign (const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& rhs)
{ }
}; // TpetraMultiVector class (specialization GO=int, NO=EpetraNode)
#endif // TpetraMultiVector class (specialization GO=long long, NO=EpetraNode)
#endif // HAVE_XPETRA_EPETRA
} // Xpetra namespace
// Following header file inculsion is needed for the dynamic_cast to TpetraVector in elementWiseMultiply (because we cannot dynamic_cast if target is not a complete type)
// It is included here to avoid circular dependency between Vector and MultiVector
// TODO: there is certainly a more elegant solution...
#include "Xpetra_TpetraVector.hpp"
namespace Xpetra {
template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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) {
XPETRA_MONITOR("TpetraMultiVector::elementWiseMultiply");
// XPETRA_DYNAMIC_CAST won't take TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
// as an argument, hence the following typedef.
typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> tpv;
XPETRA_DYNAMIC_CAST(const tpv, A, tA, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
XPETRA_DYNAMIC_CAST(const TpetraMultiVector, B, tB, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
vec_->elementWiseMultiply(scalarAB, *tA.getTpetra_Vector(), *tB.getTpetra_MultiVector(), scalarThis);
}
} // Xpetra namespace
#define XPETRA_TPETRAMULTIVECTOR_SHORT
#endif // XPETRA_TPETRAMULTIVECTOR_HPP
|