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 RCP<Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> > toTpetra(Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
00069
00070 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00071 RCP<Tpetra::Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> > toTpetra(const Vector<Scalar,LocalOrdinal, GlobalOrdinal, Node> &);
00072
00073 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00074 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
00075
00076 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00077 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec);
00078
00079
00080
00081
00082 template <class Scalar = Vector<>::scalar_type,
00083 class LocalOrdinal = typename Vector<Scalar>::local_ordinal_type,
00084 class GlobalOrdinal = typename Vector<Scalar, LocalOrdinal>::global_ordinal_type,
00085 class Node = typename Vector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
00086 class TpetraVector
00087 : public virtual Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>,
00088 public TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
00089 {
00090 #undef XPETRA_TPETRAMULTIVECTOR_SHORT
00091 #include "Xpetra_UseShortNames.hpp"
00092
00093 public:
00094
00095 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dot;
00096 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::norm1;
00097 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::norm2;
00098 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::normInf;
00099 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::normWeighted;
00100 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::meanValue;
00101 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceGlobalValue;
00102 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::sumIntoGlobalValue;
00103 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceLocalValue;
00104 using TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::sumIntoLocalValue;
00105
00107
00108
00110 TpetraVector(const Teuchos::RCP<const Map> &map, bool zeroOut=true)
00111 : TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map,1,zeroOut) { }
00112
00114 TpetraVector(const Teuchos::RCP<const Map> &map, const Teuchos::ArrayView< const Scalar > &A)
00115 : TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map,A,map->getNodeNumElements(),1) { }
00116
00118 virtual ~TpetraVector() { }
00119
00121
00123
00124
00126 void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::replaceGlobalValue"); getTpetra_Vector()->replaceGlobalValue(globalRow, value); }
00127
00129 void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::sumIntoGlobalValue"); getTpetra_Vector()->sumIntoGlobalValue(globalRow, value); }
00130
00132 void replaceLocalValue(LocalOrdinal myRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::replaceLocalValue"); getTpetra_Vector()->replaceLocalValue(myRow, value); }
00133
00135 void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value) { XPETRA_MONITOR("TpetraVector::sumIntoLocalValue"); getTpetra_Vector()->sumIntoLocalValue(myRow, value); }
00136
00138
00140
00141
00143 typename Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const { XPETRA_MONITOR("TpetraVector::norm1"); return getTpetra_Vector()->norm1(); }
00144
00146 typename Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const { XPETRA_MONITOR("TpetraVector::norm2"); return getTpetra_Vector()->norm2(); }
00147
00149 typename Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const { XPETRA_MONITOR("TpetraVector::normInf"); return getTpetra_Vector()->normInf(); }
00150
00152 Scalar meanValue() const { XPETRA_MONITOR("TpetraVector::meanValue"); return getTpetra_Vector()->meanValue(); }
00153
00155
00157
00158
00160 std::string description() const { XPETRA_MONITOR("TpetraVector::description"); return getTpetra_Vector()->description(); }
00161
00163 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraVector::describe"); getTpetra_Vector()->describe(out, verbLevel); }
00164
00166
00168 Scalar dot(const Vector &a) const { XPETRA_MONITOR("TpetraVector::dot"); return getTpetra_Vector()->dot(*toTpetra(a)); }
00169
00171 typename Teuchos::ScalarTraits< Scalar >::magnitudeType normWeighted(const Vector &weights) const { XPETRA_MONITOR("TpetraVector::normWeighted"); return getTpetra_Vector()->normWeighted(*toTpetra(weights)); }
00172
00173
00175
00176
00178 TpetraVector(const Teuchos::RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &vec) : TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (vec) { }
00179
00181 RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > getTpetra_Vector() const { return this->TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getTpetra_MultiVector()->getVectorNonConst(0); }
00182
00184
00185 };
00186
00187
00188 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00189 RCP<Tpetra::Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> > toTpetra(Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
00190 typedef TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass;
00191 XPETRA_DYNAMIC_CAST( TpetraVectorClass, x, tX, "toTpetra");
00192 return tX.getTpetra_Vector();
00193 }
00194
00195 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00196 RCP<Tpetra::Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> > toTpetra(const Vector< Scalar,LocalOrdinal, GlobalOrdinal, Node> &x) {
00197 typedef TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass;
00198 XPETRA_DYNAMIC_CAST(const TpetraVectorClass, x, tX, "toTpetra");
00199 return tX.getTpetra_Vector();
00200 }
00201
00202 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00203 RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
00204 if (!vec.is_null())
00205 return rcp(new TpetraVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(vec));
00206
00207 return Teuchos::null;
00208 }
00209
00210 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00211 RCP<const Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP<const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec) {
00212
00213 return toXpetra(Teuchos::rcp_const_cast<Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node > > (vec));
00214 }
00215
00216
00217
00218 }
00219
00220 #define XPETRA_TPETRAVECTOR_SHORT
00221 #endif // XPETRA_TPETRAVECTOR_HPP