Go to the documentation of this file.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_EPETRAVECTOR_HPP
00047 #define XPETRA_EPETRAVECTOR_HPP
00048
00049
00050
00051 #include "Xpetra_EpetraConfigDefs.hpp"
00052
00053 #include "Xpetra_Vector.hpp"
00054 #include "Xpetra_EpetraMultiVector.hpp"
00055 #include "Xpetra_EpetraMap.hpp"
00056 #include "Xpetra_Utils.hpp"
00057 #include "Xpetra_EpetraImport.hpp"
00058 #include "Xpetra_EpetraExport.hpp"
00059 #include "Xpetra_Exceptions.hpp"
00060
00061 #include <Epetra_Vector.h>
00062
00063 namespace Xpetra {
00064
00065
00066 template<class GlobalOrdinal>
00067 Epetra_Vector & toEpetra(Vector<double, int, GlobalOrdinal> &);
00068
00069 template<class GlobalOrdinal>
00070 const Epetra_Vector & toEpetra(const Vector<double, int, GlobalOrdinal> &);
00071
00072
00073 template<class EpetraGlobalOrdinal>
00074 class EpetraVectorT
00075 : public virtual Vector<double,int,EpetraGlobalOrdinal>, public EpetraMultiVectorT<EpetraGlobalOrdinal>
00076 {
00077
00078 typedef double Scalar;
00079 typedef int LocalOrdinal;
00080 typedef EpetraGlobalOrdinal GlobalOrdinal;
00081 typedef KokkosClassic::DefaultNode::DefaultNodeType Node;
00082
00083 public:
00084
00085 using EpetraMultiVectorT<GlobalOrdinal>::dot;
00086 using EpetraMultiVectorT<GlobalOrdinal>::norm1;
00087 using EpetraMultiVectorT<GlobalOrdinal>::norm2;
00088 using EpetraMultiVectorT<GlobalOrdinal>::normInf;
00089 using EpetraMultiVectorT<GlobalOrdinal>::normWeighted;
00090 using EpetraMultiVectorT<GlobalOrdinal>::meanValue;
00091 using EpetraMultiVectorT<GlobalOrdinal>::replaceGlobalValue;
00092 using EpetraMultiVectorT<GlobalOrdinal>::sumIntoGlobalValue;
00093 using EpetraMultiVectorT<GlobalOrdinal>::replaceLocalValue;
00094 using EpetraMultiVectorT<GlobalOrdinal>::sumIntoLocalValue;
00095
00097
00098
00100 explicit EpetraVectorT(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, bool zeroOut=true);
00101
00102
00104
00105
00107 virtual ~EpetraVectorT() { }
00108
00110
00112
00113
00115 void replaceGlobalValue(GlobalOrdinal globalRow, const Scalar &value);
00116
00118 void sumIntoGlobalValue(GlobalOrdinal globalRow, const Scalar &value);
00119
00121 void replaceLocalValue(LocalOrdinal myRow, const Scalar &value);
00122
00124 void sumIntoLocalValue(LocalOrdinal myRow, const Scalar &value);
00125
00127
00129
00130
00132 Scalar dot(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const;
00133
00135 Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const;
00136
00138 Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const;
00139
00141 Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const;
00142
00144 Teuchos::ScalarTraits< Scalar >::magnitudeType normWeighted(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &weights) const;
00145
00147 Scalar meanValue() const;
00148
00150
00152
00153
00155 std::string description() const;
00156
00158 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const;
00159
00161
00163
00164
00166 EpetraVectorT(const Teuchos::RCP<Epetra_Vector> &vec) : EpetraMultiVectorT<GlobalOrdinal>(vec) { }
00167
00169 Epetra_Vector * getEpetra_Vector() const { return (*this->EpetraMultiVectorT<GlobalOrdinal>::getEpetra_MultiVector())(0); }
00170
00171
00175 EpetraVectorT(const RCP<Epetra_MultiVector> &mv, size_t j);
00176
00178
00179 private:
00180
00181
00182 const RCP<const Epetra_MultiVector> internalRefToBaseMV_;
00183
00184 };
00185
00186 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
00187 typedef EpetraVectorT<int> EpetraVector;
00188 #endif
00189
00190 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
00191 typedef EpetraVectorT<long long> EpetraVector64;
00192 #endif
00193
00194 }
00195
00196 #endif // XPETRA_EPETRAVECTOR_HPP