Zoltan2
Zoltan2_XpetraMultiVectorAdapter.hpp
Go to the documentation of this file.
00001 // @HEADER
00002 //
00003 // ***********************************************************************
00004 //
00005 //   Zoltan2: A package of combinatorial algorithms for scientific computing
00006 //                  Copyright 2012 Sandia Corporation
00007 //
00008 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00009 // the U.S. Government retains certain rights in this software.
00010 //
00011 // Redistribution and use in source and binary forms, with or without
00012 // modification, are permitted provided that the following conditions are
00013 // met:
00014 //
00015 // 1. Redistributions of source code must retain the above copyright
00016 // notice, this list of conditions and the following disclaimer.
00017 //
00018 // 2. Redistributions in binary form must reproduce the above copyright
00019 // notice, this list of conditions and the following disclaimer in the
00020 // documentation and/or other materials provided with the distribution.
00021 //
00022 // 3. Neither the name of the Corporation nor the names of the
00023 // contributors may be used to endorse or promote products derived from
00024 // this software without specific prior written permission.
00025 //
00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00037 //
00038 // Questions? Contact Karen Devine      (kddevin@sandia.gov)
00039 //                    Erik Boman        (egboman@sandia.gov)
00040 //                    Siva Rajamanickam (srajama@sandia.gov)
00041 //
00042 // ***********************************************************************
00043 //
00044 // @HEADER
00045 
00050 #ifndef _ZOLTAN2_XPETRAMULTIVECTORADAPTER_HPP_
00051 #define _ZOLTAN2_XPETRAMULTIVECTORADAPTER_HPP_
00052 
00053 #include <Zoltan2_XpetraTraits.hpp>
00054 #include <Zoltan2_VectorAdapter.hpp>
00055 #include <Zoltan2_StridedData.hpp>
00056 #include <Zoltan2_Util.hpp>
00057 
00058 #include <Xpetra_EpetraMultiVector.hpp>
00059 #include <Xpetra_TpetraMultiVector.hpp>
00060 
00061 namespace Zoltan2 {
00062 
00080 template <typename User>
00081   class XpetraMultiVectorAdapter : public VectorAdapter<User> {
00082 public:
00083 
00084 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00085   typedef typename InputTraits<User>::scalar_t    scalar_t;
00086   typedef typename InputTraits<User>::lno_t    lno_t;
00087   typedef typename InputTraits<User>::gno_t    gno_t;
00088   typedef typename InputTraits<User>::zgid_t    zgid_t;
00089   typedef typename InputTraits<User>::part_t   part_t;
00090   typedef typename InputTraits<User>::node_t   node_t;
00091   typedef VectorAdapter<User>       base_adapter_t;
00092   typedef User user_t;
00093   typedef User userCoord_t;
00094 
00095   typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
00096   typedef Xpetra::TpetraMultiVector<
00097     scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
00098   typedef Xpetra::EpetraMultiVector xe_mvector_t;
00099 #endif
00100 
00103   ~XpetraMultiVectorAdapter() { }
00104 
00120   XpetraMultiVectorAdapter(const RCP<const User> &invector,
00121     std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides);
00122 
00128   XpetraMultiVectorAdapter(const RCP<const User> &invector);
00129 
00130 
00132   // The Adapter interface.
00134 
00135   size_t getLocalNumIDs() const { return vector_->getLocalLength();}
00136 
00137   void getIDsView(const zgid_t *&ids) const
00138   { 
00139     ids = map_->getNodeElementList().getRawPtr();
00140   }
00141 
00142   int getNumWeightsPerID() const { return numWeights_;}
00143 
00144   void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
00145   {
00146     env_->localInputAssertion(__FILE__, __LINE__, "invalid weight index",
00147       idx >= 0 && idx < numWeights_, BASIC_ASSERTION);
00148     size_t length;
00149     weights_[idx].getStridedList(length, weights, stride);
00150   }
00151 
00153   // The VectorAdapter interface.
00155 
00156   int getNumEntriesPerID() const {return vector_->getNumVectors();}
00157 
00158   void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const;
00159 
00160   template <typename Adapter>
00161     void applyPartitioningSolution(const User &in, User *&out,
00162          const PartitioningSolution<Adapter> &solution) const;
00163 
00164 private:
00165 
00166   RCP<const User> invector_;
00167   RCP<const x_mvector_t> vector_;
00168   RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
00169   RCP<Environment> env_;    // for error messages, etc.
00170   lno_t base_;
00171 
00172   int numWeights_;
00173   ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
00174 };
00175 
00177 // Definitions
00179 
00180 template <typename User>
00181   XpetraMultiVectorAdapter<User>::XpetraMultiVectorAdapter(
00182     const RCP<const User> &invector,
00183     std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides):
00184       invector_(invector), vector_(), map_(), 
00185       env_(rcp(new Environment)), base_(),
00186       numWeights_(weights.size()), weights_(weights.size())
00187 {
00188   typedef StridedData<lno_t, scalar_t> input_t;
00189 
00190   vector_ = XpetraTraits<User>::convertToXpetra(invector);
00191   map_ = vector_->getMap();
00192   base_ = map_->getIndexBase();
00193 
00194   size_t length = vector_->getLocalLength();
00195 
00196   if (length > 0 && numWeights_ > 0){
00197     int stride = 1;
00198     for (int w=0; w < numWeights_; w++){
00199       if (weightStrides.size())
00200         stride = weightStrides[w];
00201       ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride*length, false); 
00202       weights_[w] = input_t(wgtV, stride);
00203     }
00204   }
00205 }
00206 
00207 
00208 template <typename User>
00209   XpetraMultiVectorAdapter<User>::XpetraMultiVectorAdapter(
00210     const RCP<const User> &invector):
00211       invector_(invector), vector_(), map_(), 
00212       env_(rcp(new Environment)), base_(),
00213       numWeights_(0), weights_()
00214 {
00215   vector_ = XpetraTraits<User>::convertToXpetra(invector);
00216   map_ = vector_->getMap();
00217   base_ = map_->getIndexBase();
00218 }
00219 
00220 template <typename User>
00221   void XpetraMultiVectorAdapter<User>::getEntriesView(
00222     const scalar_t *&elements, int &stride, int idx) const
00223 {
00224   size_t vecsize;
00225   stride = 1;
00226   elements = NULL;
00227   if (map_->lib() == Xpetra::UseTpetra){
00228     const xt_mvector_t *tvector = 
00229       dynamic_cast<const xt_mvector_t *>(vector_.get());
00230      
00231     vecsize = tvector->getLocalLength();
00232     if (vecsize > 0){
00233       ArrayRCP<const scalar_t> data = tvector->getData(idx);
00234       elements = data.get();
00235     }
00236   }
00237   else if (map_->lib() == Xpetra::UseEpetra){
00238     const xe_mvector_t *evector = 
00239       dynamic_cast<const xe_mvector_t *>(vector_.get());
00240       
00241     vecsize = evector->getLocalLength();
00242     if (vecsize > 0){
00243       ArrayRCP<const double> data = evector->getData(idx);
00244 
00245       // Cast so this will compile when scalar_t is not double,
00246       // a case when this code should never execute.
00247       elements = reinterpret_cast<const scalar_t *>(data.get());
00248     }
00249   }
00250   else{
00251     throw std::logic_error("invalid underlying lib");
00252   }
00253 }
00254 
00255 template <typename User>
00256   template <typename Adapter>
00257     void XpetraMultiVectorAdapter<User>::applyPartitioningSolution(
00258       const User &in, User *&out, 
00259       const PartitioningSolution<Adapter> &solution) const
00260 {
00261   size_t len = solution.getLocalNumberOfIds();
00262   const zgid_t *gids = solution.getIdList();
00263   const part_t *parts = solution.getPartList();
00264   ArrayRCP<zgid_t> gidList = arcp(const_cast<zgid_t *>(gids), 0, len, false);
00265   ArrayRCP<part_t> partList = arcp(const_cast<part_t *>(parts), 0, len, 
00266     false);
00267   ArrayRCP<lno_t> dummyIn;
00268   ArrayRCP<zgid_t> importList;
00269   ArrayRCP<lno_t> dummyOut;
00270   size_t numNewRows;
00271 
00272   const RCP<const Comm<int> > comm = map_->getComm();
00273 
00274   try{
00275     // Get an import list
00276     numNewRows = solution.convertSolutionToImportList(
00277       0, dummyIn, importList, dummyOut);
00278   }
00279   Z2_FORWARD_EXCEPTIONS;
00280 
00281   RCP<const User> inPtr = rcp(&in, false);
00282 
00283   RCP<const User> outPtr = XpetraTraits<User>::doMigration(
00284    inPtr, numNewRows, importList.get());
00285 
00286   out = const_cast<User *>(outPtr.get());
00287   outPtr.release();
00288 }
00289   
00290 }  //namespace Zoltan2
00291   
00292 #endif