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_TPETRAVECTOR_HPP
00047 #define XPETRA_TPETRAVECTOR_HPP
00048
00049
00050
00051 #include "Xpetra_TpetraConfigDefs.hpp"
00052
00053 #include "Xpetra_Vector.hpp"
00054 #include "Xpetra_MultiVector.hpp"
00055 #include "Xpetra_TpetraMultiVector.hpp"
00056
00057 #include "Xpetra_TpetraMap.hpp"
00058 #include "Xpetra_Utils.hpp"
00059 #include "Xpetra_TpetraImport.hpp"
00060 #include "Xpetra_TpetraExport.hpp"
00061
00062 #include "Tpetra_Vector.hpp"
00063
00064 namespace Xpetra {
00065
00066
00067 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00068 Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
00069
00070 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00071 const Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(const Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
00072
00073
00074 template <class Scalar, class LocalOrdinal = int, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
00075 class TpetraVector
00076 : public virtual Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>, public TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
00077 {
00078
00079 public:
00080
00081 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dot;
00082 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm1;
00083 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm2;
00084 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::normInf;
00085 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::normWeighted;
00086 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::meanValue;
00087 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValue;
00088 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalValue;
00089 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValue;
00090 using TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValue;
00091
00093
00094
00096 TpetraVector(const Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> > &map, bool zeroOut=true)
00097 : TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(map,1,zeroOut) { }
00098
00100 TpetraVector(const Teuchos::RCP< const Map<LocalOrdinal,GlobalOrdinal,Node> > &map, const Teuchos::ArrayView< const Scalar > &A)
00101 : TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(map,A,map->getNodeNumElements(),1) { }
00102
00104 virtual ~TpetraVector() { }
00105
00107
00109
00110
00112 void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::replaceGlobalValue"); getTpetra_Vector()->replaceGlobalValue(globalRow, value); }
00113
00115 void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::sumIntoGlobalValue"); getTpetra_Vector()->sumIntoGlobalValue(globalRow, value); }
00116
00118 void replaceLocalValue(LocalOrdinal myRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::replaceLocalValue"); getTpetra_Vector()->replaceLocalValue(myRow, value); }
00119
00121 void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::sumIntoLocalValue"); getTpetra_Vector()->sumIntoLocalValue(myRow, value); }
00122
00124
00126
00127
00129 Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const { XPETRA_MONITOR("TpetraVector::dot"); return getTpetra_Vector()->dot(toTpetra(a)); }
00130
00132 typename Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const { XPETRA_MONITOR("TpetraVector::norm1"); return getTpetra_Vector()->norm1(); }
00133
00135 typename Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const { XPETRA_MONITOR("TpetraVector::norm2"); return getTpetra_Vector()->norm2(); }
00136
00138 typename Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const { XPETRA_MONITOR("TpetraVector::normInf"); return getTpetra_Vector()->normInf(); }
00139
00141 typename Teuchos::ScalarTraits< Scalar >::magnitudeType normWeighted(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &weights) const { XPETRA_MONITOR("TpetraVector::normWeighted"); return getTpetra_Vector()->normWeighted(toTpetra(weights)); }
00142
00144 Scalar meanValue() const { XPETRA_MONITOR("TpetraVector::meanValue"); return getTpetra_Vector()->meanValue(); }
00145
00147
00149
00150
00152 std::string description() const { XPETRA_MONITOR("TpetraVector::description"); return getTpetra_Vector()->description(); }
00153
00155 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraVector::describe"); getTpetra_Vector()->describe(out, verbLevel); }
00156
00158
00160
00161
00163 TpetraVector(const Teuchos::RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) : TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(vec) { }
00164
00166 RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_Vector() const { return this->TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getTpetra_MultiVector()->getVectorNonConst(0); }
00167
00169
00170 };
00171
00172
00173 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00174 Tpetra::Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
00175 typedef TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass;
00176 XPETRA_DYNAMIC_CAST( TpetraVectorClass, x, tX, "toTpetra");
00177 return *tX.getTpetra_Vector();
00178 }
00179
00180 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00181 const Tpetra::Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> & toTpetra(const Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
00182 typedef TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass;
00183 XPETRA_DYNAMIC_CAST(const TpetraVectorClass, x, tX, "toTpetra");
00184 return *tX.getTpetra_Vector();
00185 }
00186
00187
00188 }
00189
00190 #define XPETRA_TPETRAVECTOR_SHORT
00191 #endif // XPETRA_TPETRAVECTOR_HPP