Zoltan2 Version of the Day
Zoltan2_XpetraVectorInput.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_XPETRAVECTORINPUT_HPP_
00051 #define _ZOLTAN2_XPETRAVECTORINPUT_HPP_
00052 
00053 #include <Zoltan2_VectorInput.hpp>
00054 #include <Zoltan2_XpetraTraits.hpp>
00055 #include <Zoltan2_StridedData.hpp>
00056 #include <Zoltan2_Util.hpp>
00057 
00058 #include <Xpetra_EpetraVector.hpp>
00059 #include <Xpetra_TpetraVector.hpp>
00060 
00061 namespace Zoltan2 {
00062 
00063 
00083 template <typename User>
00084   class XpetraVectorInput : public VectorInput<User> {
00085 public:
00086 
00087 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00088   typedef typename InputTraits<User>::scalar_t    scalar_t;
00089   typedef typename InputTraits<User>::lno_t    lno_t;
00090   typedef typename InputTraits<User>::gno_t    gno_t;
00091   typedef typename InputTraits<User>::gid_t    gid_t;
00092   typedef typename InputTraits<User>::node_t   node_t;
00093   typedef VectorInput<User>       base_adapter_t;
00094   typedef User user_t;
00095 
00096   typedef Xpetra::Vector<
00097     scalar_t, lno_t, gno_t, node_t> x_vector_t;
00098   typedef Xpetra::TpetraVector<
00099     scalar_t, lno_t, gno_t, node_t> xt_vector_t;
00100   typedef Xpetra::EpetraVector xe_vector_t;
00101 #endif
00102 
00105   ~XpetraVectorInput() { }
00106 
00121   XpetraVectorInput( const RCP<const User> &invector,
00122     vector<const scalar_t *> &weights, vector<int> &weightStrides);
00123 
00127   const RCP<const x_vector_t> &getVector() const
00128   {
00129     return vector_;
00130   }
00131 
00133   // The InputAdapter interface.
00135 
00136   string inputAdapterName()const { return string("XpetraVector");}
00137 
00138   size_t getLocalNumberOfObjects() const { return getLocalLength();}
00139 
00140   int getNumberOfWeightsPerObject() const { return numWeights_;}
00141 
00142   size_t getObjectWeights(int dim, const scalar_t *&wgt, int &stride) const
00143   {
00144     return getVectorWeights(dim, wgt, stride);
00145   }
00146 
00148   // The VectorInput interface.
00150 
00151   int getNumberOfVectors() const { return 1; }
00152 
00153   int getNumberOfWeights() const {return numWeights_;}
00154 
00155   size_t getLocalLength() const {return vector_->getLocalLength();}
00156   
00157   size_t getGlobalLength() const {return vector_->getGlobalLength();}
00158 
00159   size_t getVector(const gid_t *&Ids, const scalar_t *&elements, 
00160     int &stride) const;
00161 
00162   size_t getVector(int vectorNumber, const gid_t *&Ids, 
00163     const scalar_t *&elements, int &stride) const
00164   {
00165     env_->localInputAssertion(__FILE__, __LINE__, "invalid vector",
00166       vectorNumber==0, BASIC_ASSERTION);
00167 
00168     return getVector(Ids, elements, stride);
00169   }
00170 
00171   size_t getVectorWeights(int dim, const scalar_t *&weights, int &stride) const
00172   {
00173     env_->localInputAssertion(__FILE__, __LINE__, "invalid dimension",
00174       dim >= 0 && dim < numWeights_, BASIC_ASSERTION);
00175 
00176     size_t length;
00177 
00178     weights_[dim].getStridedList(length, weights, stride);
00179 
00180     return length;
00181   }
00182 
00183   template <typename Adapter>
00184     size_t applyPartitioningSolution(const User &in, User *&out,
00185          const PartitioningSolution<Adapter> &solution) const;
00186 
00187 private:
00188 
00189   RCP<const User> invector_;
00190   RCP<const x_vector_t> vector_;
00191   RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
00192   RCP<Environment> env_;
00193   lno_t base_;
00194 
00195   int numWeights_;
00196   ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
00197 };
00198 
00200 // Definitions
00202   
00203 template <typename User>
00204   XpetraVectorInput<User>::XpetraVectorInput(
00205     const RCP<const User> &invector, 
00206     vector<const scalar_t *> &weights, vector<int> &weightStrides):
00207       invector_(invector), vector_(), map_(),
00208       env_(rcp(new Environment)), base_(),
00209       numWeights_(weights.size()), weights_()
00210 {
00211   typedef StridedData<lno_t, scalar_t> input_t;
00212 
00213   vector_ = XpetraTraits<User>::convertToXpetra(invector);
00214   map_ = vector_->getMap();
00215   base_ = map_->getIndexBase();
00216 
00217   size_t length = vector_->getLocalLength();
00218 
00219   if (numWeights_ > 0)
00220     weights_ = arcp(new input_t [numWeights_], 0, numWeights_, true);
00221 
00222   if (length > 0 && numWeights_ > 0){
00223     int stride = 1;
00224     for (int w=0; w < numWeights_; w++){
00225       if (weightStrides.size())
00226         stride = weightStrides[w];
00227       ArrayRCP<const scalar_t> wtArray(weights[w], 0, stride*length, false);
00228       weights_[w] = input_t(wtArray, stride);
00229     }
00230   }
00231 }
00232 
00233 template <typename User>
00234   size_t XpetraVectorInput<User>::getVector(const gid_t *&Ids, 
00235     const scalar_t *&elements, int &stride) const
00236 {
00237   stride = 1;
00238   elements = NULL;
00239   const x_vector_t *vec =  vector_.get();
00240 
00241   if (map_->lib() == Xpetra::UseTpetra){
00242     const xt_vector_t *tvector = dynamic_cast<const xt_vector_t *>(vec);
00243 
00244     if (tvector->getLocalLength() > 0){
00245       // getData hangs if vector length is 0
00246       ArrayRCP<const scalar_t> data = tvector->getData(0);
00247       elements = data.get();
00248     }
00249   }
00250   else if (map_->lib() == Xpetra::UseEpetra){
00251     const xe_vector_t *evector = dynamic_cast<const xe_vector_t *>(vec);
00252       
00253     if (evector->getLocalLength() > 0){
00254       // getData hangs if vector length is 0
00255       ArrayRCP<const double> data = evector->getData(0);
00256 
00257       // Cast so this will compile when scalar_t is not double,
00258       // a case when this code should never execute.
00259       elements = reinterpret_cast<const scalar_t *>(data.get());
00260     }
00261   }
00262   else{
00263     throw logic_error("invalid underlying lib");
00264   }
00265 
00266   ArrayView<const gid_t> gids = map_->getNodeElementList();
00267   Ids = gids.getRawPtr();
00268 
00269   return getLocalLength();
00270 }
00271 
00272 template <typename User>
00273   template <typename Adapter>
00274     size_t XpetraVectorInput<User>::applyPartitioningSolution(
00275       const User &in, User *&out, 
00276       const PartitioningSolution<Adapter> &solution) const
00277 { 
00278   // Get an import list
00279 
00280   size_t len = solution.getLocalNumberOfIds();
00281   const gid_t *gids = solution.getIdList();
00282   const partId_t *parts = solution.getPartList();
00283   ArrayRCP<gid_t> gidList = arcp(const_cast<gid_t *>(gids), 0, len, false);
00284   ArrayRCP<partId_t> partList = arcp(const_cast<partId_t *>(parts), 0, len, 
00285     false);
00286 
00287   ArrayRCP<lno_t> dummyIn;
00288   ArrayRCP<gid_t> importList;
00289   ArrayRCP<lno_t> dummyOut;
00290   size_t numNewRows;
00291 
00292   const RCP<const Comm<int> > comm = map_->getComm(); 
00293 
00294   try{
00295     numNewRows = solution.convertSolutionToImportList(
00296       0, dummyIn, importList, dummyOut);
00297   }
00298   Z2_FORWARD_EXCEPTIONS;
00299 
00300   RCP<const User> inPtr = rcp(&in, false);
00301   lno_t localNumElts = numNewRows;
00302 
00303   RCP<const User> outPtr = XpetraTraits<User>::doMigration(
00304    inPtr, localNumElts, importList.get());
00305 
00306   out = const_cast<User *>(outPtr.get());
00307   outPtr.release();
00308   return numNewRows;
00309 }
00310 
00311 }  //namespace Zoltan2
00312   
00313 #endif