PlayaMPIContainerComm.hpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                 Playa: Programmable Linear Algebra
00005 //                 Copyright 2012 Sandia Corporation
00006 // 
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Kevin Long (kevin.long@ttu.edu)
00038 // 
00039 
00040 /* @HEADER@ */
00041 
00042 #ifndef PLAYA_MPICONTAINERCOMM_H
00043 #define PLAYA_MPICONTAINERCOMM_H
00044 
00045 /*! \file PlayaMPIContainerComm.hpp
00046   \brief Object representation of an MPI communicator for templated containers
00047 */
00048 
00049 #include "PlayaMPIComm.hpp"
00050 #include "PlayaMPITraits.hpp"
00051 
00052 namespace Playa
00053 {
00054 using Teuchos::Array;
00055   /** \ingroup MPI
00056    * \brief Object representation of an MPI communicator for templated containers
00057    * \note Template specialization exists for <tt>std::string</tt>.
00058    * @author Kevin Long
00059    */
00060 
00061   template <class T> class MPIContainerComm
00062   {
00063   public:
00064 
00065     //! Broadcast a single object 
00066     static void bcast(T& x, int src, const MPIComm& comm);
00067 
00068     //! Broadcast an array of objects
00069     static void bcast(Array<T>& x, int src, const MPIComm& comm);
00070 
00071     //! Broadcast an array of arrays 
00072     static void bcast(Array<Array<T> >& x,
00073                       int src, const MPIComm& comm);
00074 
00075     //! Gather to all processors
00076     static void allGather(const T& outgoing,
00077                           Array<T>& incoming,
00078                           const MPIComm& comm);
00079 
00080     //! All-to-all scatter/gather for an array of objects
00081     static void allToAll(const Array<T>& outgoing,
00082                          Array<Array<T> >& incoming,
00083                          const MPIComm& comm);
00084 
00085     //! All-to-all scatter/gather for an array of arrays
00086     static void allToAll(const Array<Array<T> >& outgoing,
00087                          Array<Array<T> >& incoming,
00088                          const MPIComm& comm);
00089 
00090     /** Gatherv: gather arrays of data to the root processor */
00091     static void gatherv(const Array<T>& outgoing,
00092                         Array<Array<T> >& incoming,
00093                         int rootRank,
00094                         const MPIComm& comm);
00095 
00096     //! Sum local values from all processors with rank < myRank
00097     static void accumulate(const T& localValue, Array<T>& sums, T& total,
00098                            const MPIComm& comm);
00099 
00100   private:
00101     //! Build a 1D array and an offset list from a 2D array
00102     static void getBigArray(const Array<Array<T> >& x,
00103                             Array<T>& bigArray,
00104                             Array<int>& offsets);
00105 
00106     //! Reassemble a 2D array from a 1D array and an offset table
00107     static void getSmallArrays(const Array<T>& bigArray,
00108                                const Array<int>& offsets,
00109                                Array<Array<T> >& x);
00110 
00111 
00112   };
00113 
00114 
00115 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00116   /** \ingroup MPI
00117    * Specialiaztion of MPIContainerComm<T> to std::string
00118    */
00119   template <> class MPIContainerComm<std::string>
00120   {
00121   public:
00122     static void bcast(std::string& x, int src, const MPIComm& comm);
00123 
00124     /** bcast an array of objects */
00125     static void bcast(Array<std::string>& x, int src, const MPIComm& comm);
00126 
00127     /** bcast an array of arrays  */
00128     static void bcast(Array<Array<std::string> >& x,
00129                       int src, const MPIComm& comm);
00130 
00131     /** AllGather: each process sends a single object to all other procs */
00132     static void allGather(const std::string& outgoing,
00133                           Array<std::string>& incoming,
00134                           const MPIComm& comm);
00135 
00136     /** Gatherv: gather arrays of strings to the root processor */
00137     static void gatherv(const Array<std::string>& outgoing,
00138                         Array<Array<std::string> >& incoming,
00139                         int rootRank,
00140                         const MPIComm& comm);
00141 
00142     /** get a single big array of characters from an array of strings,
00143      * packing the structural description into the header of the packed
00144      * array as follows:
00145      * \code
00146      * [numStrings, offset0, offset1, ..., offsetN, char data]
00147      * \endcode 
00148      */
00149     static void pack(const Array<std::string>& x,
00150                      Array<char>& packed);
00151 
00152     /** recover an array of strings from a single big array and
00153      * and offset table */
00154     static void unpack(const Array<char>& packed,
00155                        Array<std::string>& x);
00156   private:
00157     /** get a single big array of characters from an array of strings */
00158     static void getBigArray(const Array<std::string>& x,
00159                             Array<char>& bigArray,
00160                             Array<int>& offsets);
00161 
00162     /** recover an array of strings from a single big array and
00163      * and offset table */
00164     static void getStrings(const Array<char>& bigArray,
00165                            const Array<int>& offsets,
00166                            Array<std::string>& x);
00167   };
00168 
00169 #endif // DOXYGEN_SHOULD_SKIP_THIS
00170 
00171   /* --------- generic functions for primitives ------------------- */
00172 
00173   template <class T> inline void MPIContainerComm<T>::bcast(T& x, int src,
00174                                                             const MPIComm& comm)
00175   {
00176     comm.bcast((void*)&x, 1, MPITraits<T>::type(), src);
00177   }
00178 
00179 
00180   /* ----------- generic functions for arrays of primitives ----------- */
00181 
00182   template <class T>
00183   inline void MPIContainerComm<T>::bcast(Array<T>& x, int src, const MPIComm& comm)
00184   {
00185     int len = x.length();
00186     MPIContainerComm<int>::bcast(len, src, comm);
00187 
00188     if (comm.getRank() != src)
00189       {
00190         x.resize(len);
00191       }
00192     if (len==0) return;
00193 
00194     /* then broadcast the contents */
00195     comm.bcast((void*) &(x[0]), (int) len,
00196                MPITraits<T>::type(), src);
00197   }
00198 
00199 
00200 
00201   /* ---------- generic function for arrays of arrays ----------- */
00202 
00203   template <class T>
00204   inline void MPIContainerComm<T>::bcast(Array<Array<T> >& x, int src, const MPIComm& comm)
00205   {
00206     Array<T> bigArray;
00207     Array<int> offsets;
00208 
00209     if (src==comm.getRank())
00210       {
00211         getBigArray(x, bigArray, offsets);
00212       }
00213 
00214     bcast(bigArray, src, comm);
00215     MPIContainerComm<int>::bcast(offsets, src, comm);
00216 
00217     if (src != comm.getRank())
00218       {
00219         getSmallArrays(bigArray, offsets, x);
00220       }
00221   }
00222 
00223   /* ---------- generic gather and scatter ------------------------ */
00224 
00225   template <class T> inline
00226   void MPIContainerComm<T>::allToAll(const Array<T>& outgoing,
00227                                      Array<Array<T> >& incoming,
00228                                      const MPIComm& comm)
00229   {
00230     int numProcs = comm.getNProc();
00231 
00232     // catch degenerate case
00233     if (numProcs==1)
00234       {
00235         incoming.resize(1);
00236         incoming[0] = outgoing;
00237         return;
00238       }
00239 
00240     Array<T> sb(numProcs * outgoing.length());
00241     Array<T> rb(numProcs * outgoing.length());
00242 
00243     T* sendBuf = new T[numProcs * outgoing.length()];
00244     TEUCHOS_TEST_FOR_EXCEPTION(sendBuf==0, 
00245                        std::runtime_error, "Comm::allToAll failed to allocate sendBuf");
00246 
00247     T* recvBuf = new T[numProcs * outgoing.length()];
00248     TEUCHOS_TEST_FOR_EXCEPTION(recvBuf==0, 
00249                        std::runtime_error, "Comm::allToAll failed to allocate recvBuf");
00250 
00251     int i;
00252     for (i=0; i<numProcs; i++)
00253       {
00254         for (int j=0; j<outgoing.length(); j++)
00255           {
00256             sendBuf[i*outgoing.length() + j] = outgoing[j];
00257           }
00258       }
00259 
00260 
00261 
00262     comm.allToAll(sendBuf, outgoing.length(), MPITraits<T>::type(),
00263                   recvBuf, outgoing.length(), MPITraits<T>::type());
00264 
00265     incoming.resize(numProcs);
00266 
00267     for (i=0; i<numProcs; i++)
00268       {
00269         incoming[i].resize(outgoing.length());
00270         for (int j=0; j<outgoing.length(); j++)
00271           {
00272             incoming[i][j] = recvBuf[i*outgoing.length() + j];
00273           }
00274       }
00275 
00276     delete [] sendBuf;
00277     delete [] recvBuf;
00278   }
00279 
00280   template <class T> inline
00281   void MPIContainerComm<T>::allToAll(const Array<Array<T> >& outgoing,
00282                                      Array<Array<T> >& incoming, const MPIComm& comm)
00283   {
00284     int numProcs = comm.getNProc();
00285 
00286     // catch degenerate case
00287     if (numProcs==1)
00288       {
00289         incoming = outgoing;
00290         return;
00291       }
00292 
00293     int* sendMesgLength = new int[numProcs];
00294     TEUCHOS_TEST_FOR_EXCEPTION(sendMesgLength==0, 
00295                        std::runtime_error, "failed to allocate sendMesgLength");
00296     int* recvMesgLength = new int[numProcs];
00297     TEUCHOS_TEST_FOR_EXCEPTION(recvMesgLength==0, 
00298                        std::runtime_error, "failed to allocate recvMesgLength");
00299 
00300     int p = 0;
00301     for (p=0; p<numProcs; p++)
00302       {
00303         sendMesgLength[p] = outgoing[p].length();
00304       }
00305     
00306     comm.allToAll(sendMesgLength, 1, MPIDataType::intType(),
00307                   recvMesgLength, 1, MPIDataType::intType());
00308 
00309 
00310     int totalSendLength = 0;
00311     int totalRecvLength = 0;
00312     for (p=0; p<numProcs; p++)
00313       {
00314         totalSendLength += sendMesgLength[p];
00315         totalRecvLength += recvMesgLength[p];
00316       }
00317 
00318     T* sendBuf = new T[totalSendLength];
00319     TEUCHOS_TEST_FOR_EXCEPTION(sendBuf==0, 
00320                        std::runtime_error, "failed to allocate sendBuf");
00321     T* recvBuf = new T[totalRecvLength];
00322     TEUCHOS_TEST_FOR_EXCEPTION(recvBuf==0, 
00323                        std::runtime_error, "failed to allocate recvBuf");
00324 
00325     int* sendDisp = new int[numProcs];
00326     TEUCHOS_TEST_FOR_EXCEPTION(sendDisp==0, 
00327                        std::runtime_error, "failed to allocate sendDisp");
00328     int* recvDisp = new int[numProcs];
00329     TEUCHOS_TEST_FOR_EXCEPTION(recvDisp==0, 
00330                        std::runtime_error, "failed to allocate recvDisp");
00331 
00332     int count = 0;
00333     sendDisp[0] = 0;
00334     recvDisp[0] = 0;
00335 
00336     for (p=0; p<numProcs; p++)
00337       {
00338         for (int i=0; i<outgoing[p].length(); i++)
00339           {
00340             sendBuf[count] = outgoing[p][i];
00341             count++;
00342           }
00343         if (p>0)
00344           {
00345             sendDisp[p] = sendDisp[p-1] + sendMesgLength[p-1];
00346             recvDisp[p] = recvDisp[p-1] + recvMesgLength[p-1];
00347           }
00348       }
00349 
00350     comm.allToAllv(sendBuf, sendMesgLength,
00351                    sendDisp, MPITraits<T>::type(),
00352                    recvBuf, recvMesgLength,
00353                    recvDisp, MPITraits<T>::type());
00354 
00355     incoming.resize(numProcs);
00356     for (p=0; p<numProcs; p++)
00357       {
00358         incoming[p].resize(recvMesgLength[p]);
00359         for (int i=0; i<recvMesgLength[p]; i++)
00360           {
00361             incoming[p][i] = recvBuf[recvDisp[p] + i];
00362           }
00363       }
00364 
00365     delete [] sendBuf;
00366     delete [] sendMesgLength;
00367     delete [] sendDisp;
00368     delete [] recvBuf;
00369     delete [] recvMesgLength;
00370     delete [] recvDisp;
00371   }
00372 
00373   template <class T> inline
00374   void MPIContainerComm<T>::allGather(const T& outgoing, Array<T>& incoming,
00375                                       const MPIComm& comm)
00376   {
00377     int nProc = comm.getNProc();
00378     incoming.resize(nProc);
00379 
00380     if (nProc==1)
00381       {
00382         incoming[0] = outgoing;
00383       }
00384     else
00385       {
00386         comm.allGather((void*) &outgoing, 1, MPITraits<T>::type(),
00387                        (void*) &(incoming[0]), 1, MPITraits<T>::type());
00388       }
00389   }
00390 
00391   template <class T> inline
00392   void MPIContainerComm<T>::accumulate(const T& localValue, Array<T>& sums,
00393                                        T& total,
00394                                        const MPIComm& comm)
00395   {
00396     Array<T> contributions;
00397     allGather(localValue, contributions, comm);
00398     sums.resize(comm.getNProc());
00399     sums[0] = 0;
00400     total = contributions[0];
00401 
00402     for (int i=0; i<comm.getNProc()-1; i++)
00403       {
00404         total += contributions[i+1];
00405         sums[i+1] = sums[i] + contributions[i];
00406       }
00407   }
00408 
00409 
00410 
00411 
00412   template <class T> inline
00413   void MPIContainerComm<T>::getBigArray(const Array<Array<T> >& x, Array<T>& bigArray,
00414                                         Array<int>& offsets)
00415   {
00416     offsets.resize(x.length()+1);
00417     int totalLength = 0;
00418 
00419     for (int i=0; i<x.length(); i++)
00420       {
00421         offsets[i] = totalLength;
00422         totalLength += x[i].length();
00423       }
00424     offsets[x.length()] = totalLength;
00425 
00426     bigArray.resize(totalLength);
00427 
00428     for (int i=0; i<x.length(); i++)
00429       {
00430         for (int j=0; j<x[i].length(); j++)
00431           {
00432             bigArray[offsets[i]+j] = x[i][j];
00433           }
00434       }
00435   }
00436 
00437   template <class T> inline
00438   void MPIContainerComm<T>::getSmallArrays(const Array<T>& bigArray,
00439                                            const Array<int>& offsets,
00440                                            Array<Array<T> >& x)
00441   {
00442     x.resize(offsets.length()-1);
00443     for (int i=0; i<x.length(); i++)
00444       {
00445         x[i].resize(offsets[i+1]-offsets[i]);
00446         for (int j=0; j<x[i].length(); j++)
00447           {
00448             x[i][j] = bigArray[offsets[i] + j];
00449           }
00450       }
00451   }
00452 
00453 
00454 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00455 
00456   /* --------------- std::string specializations --------------------- */
00457 
00458   inline void MPIContainerComm<std::string>::bcast(std::string& x,
00459                                               int src, const MPIComm& comm)
00460   {
00461     int len = x.length();
00462     MPIContainerComm<int>::bcast(len, src, comm);
00463 
00464     x.resize(len);
00465     comm.bcast((void*)&(x[0]), len, MPITraits<char>::type(), src);
00466   }
00467 
00468 
00469   inline void MPIContainerComm<std::string>::bcast(Array<std::string>& x, int src,
00470                                               const MPIComm& comm)
00471   {
00472     /* begin by packing all the data into a big char array. This will
00473      * take a little time, but will be cheaper than multiple MPI calls */
00474     Array<char> bigArray;
00475     Array<int> offsets;
00476     if (comm.getRank()==src)
00477       {
00478         getBigArray(x, bigArray, offsets);
00479       }
00480 
00481     /* now broadcast the big array and the offsets */
00482     MPIContainerComm<char>::bcast(bigArray, src, comm);
00483     MPIContainerComm<int>::bcast(offsets, src, comm);
00484 
00485     /* finally, reassemble the array of strings */
00486     if (comm.getRank() != src)
00487       {
00488         getStrings(bigArray, offsets, x);
00489       }
00490   }
00491 
00492   inline void MPIContainerComm<std::string>::bcast(Array<Array<std::string> >& x,
00493                                               int src, const MPIComm& comm)
00494   {
00495     int len = x.length();
00496     MPIContainerComm<int>::bcast(len, src, comm);
00497 
00498     x.resize(len);
00499     for (int i=0; i<len; i++)
00500       {
00501         MPIContainerComm<std::string>::bcast(x[i], src, comm);
00502       }
00503   }
00504 
00505 
00506   inline void MPIContainerComm<std::string>::allGather(const std::string& outgoing,
00507                                                   Array<std::string>& incoming,
00508                                                   const MPIComm& comm)
00509   {
00510     int nProc = comm.getNProc();
00511 
00512     int sendCount = outgoing.length();
00513 
00514     incoming.resize(nProc);
00515 
00516     int* recvCounts = new int[nProc];
00517     int* recvDisplacements = new int[nProc];
00518 
00519     /* share lengths with all procs */
00520     comm.allGather((void*) &sendCount, 1, MPIDataType::intType(),
00521                    (void*) recvCounts, 1, MPIDataType::intType());
00522 
00523 
00524     int recvSize = 0;
00525     recvDisplacements[0] = 0;
00526     for (int i=0; i<nProc; i++)
00527       {
00528         recvSize += recvCounts[i];
00529         if (i < nProc-1)
00530           {
00531             recvDisplacements[i+1] = recvDisplacements[i]+recvCounts[i];
00532           }
00533       }
00534 
00535     char* recvBuf = new char[recvSize];
00536 
00537     comm.allGatherv((void*) outgoing.c_str(), sendCount, MPIDataType::charType(),
00538                     recvBuf, recvCounts, recvDisplacements, MPIDataType::charType());
00539 
00540     for (int j=0; j<nProc; j++)
00541       {
00542         char* start = recvBuf + recvDisplacements[j];
00543         char* tmp = new char[recvCounts[j]+1];
00544         std::memcpy(tmp, start, recvCounts[j]);
00545         tmp[recvCounts[j]] = '\0';
00546         incoming[j] = std::string(tmp);
00547         delete [] tmp;
00548       }
00549     
00550     delete [] recvCounts;
00551     delete [] recvDisplacements;
00552     delete [] recvBuf;
00553   }
00554   
00555   inline void MPIContainerComm<std::string>::gatherv(const Array<std::string>& outgoing,
00556                                                 Array<Array<std::string> >& incoming,
00557                                                 int root,
00558                                                 const MPIComm& comm)
00559   {
00560     int nProc = comm.getNProc();
00561 
00562     Array<char> packedLocalArray;
00563     pack(outgoing, packedLocalArray);
00564 
00565     int sendCount = packedLocalArray.size();
00566 
00567     /* gather the message sizes from all procs */
00568     Array<int> recvCounts(nProc);
00569     Array<int> recvDisplacements(nProc);
00570 
00571     comm.gather((void*) &sendCount, 1, MPIDataType::intType(),
00572                 (void*) &(recvCounts[0]), 1, MPIDataType::intType(), root);
00573     
00574     /* compute the displacements */
00575     int recvSize = 0;
00576     if (root == comm.getRank())
00577       {
00578         recvDisplacements[0] = 0;
00579         for (int i=0; i<nProc; i++)
00580           {
00581             recvSize += recvCounts[i];
00582             if (i < nProc-1)
00583               {
00584                 recvDisplacements[i+1] = recvDisplacements[i]+recvCounts[i];
00585               }
00586           }
00587       }
00588 
00589     /* set the size to 1 on non-root procs */
00590     Array<char> recvBuf(std::max(1,recvSize));
00591     
00592 
00593     void* sendBuf = (void*) &(packedLocalArray[0]);
00594     void* inBuf = (void*) &(recvBuf[0]);
00595     int* inCounts = &(recvCounts[0]);
00596     int* inDisps = &(recvDisplacements[0]);
00597 
00598     /* gather the packed data */
00599     comm.gatherv( sendBuf, sendCount, MPIDataType::charType(),
00600                   inBuf, inCounts, inDisps,
00601                   MPIDataType::charType(), root);
00602 
00603     /* on the root, unpack the data */
00604     if (comm.getRank()==root)
00605       {
00606         incoming.resize(nProc);
00607         for (int j=0; j<nProc; j++)
00608           {
00609             char* start = &(recvBuf[0]) + recvDisplacements[j];
00610             Array<char> tmp(recvCounts[j]+1);
00611             std::memcpy(&(tmp[0]), start, recvCounts[j]);
00612             tmp[recvCounts[j]] = '\0';
00613             unpack(tmp, incoming[j]);
00614           }
00615       }
00616                  
00617                  
00618   }
00619 
00620 
00621   inline void MPIContainerComm<std::string>::getBigArray(const Array<std::string>& x,
00622                                                     Array<char>& bigArray,
00623                                                     Array<int>& offsets)
00624   {
00625     offsets.resize(x.length()+1);
00626     int totalLength = 0;
00627 
00628     for (int i=0; i<x.length(); i++)
00629       {
00630         offsets[i] = totalLength;
00631         totalLength += x[i].length();
00632       }
00633     offsets[x.length()] = totalLength;
00634 
00635     bigArray.resize(totalLength);
00636 
00637     for (int i=0; i<x.length(); i++)
00638       {
00639         for (unsigned int j=0; j<x[i].length(); j++)
00640           {
00641             bigArray[offsets[i]+j] = x[i][j];
00642           }
00643       }
00644   }
00645 
00646   inline void MPIContainerComm<std::string>::pack(const Array<std::string>& x,
00647                                              Array<char>& bigArray)
00648   {
00649     Array<int> offsets(x.size()+1);
00650     int headerSize = (x.size()+2) * sizeof(int);
00651 
00652     int totalLength = headerSize;
00653 
00654     for (int i=0; i<x.length(); i++)
00655       {
00656         offsets[i] = totalLength;
00657         totalLength += x[i].length();
00658       }
00659     offsets[x.length()] = totalLength;
00660 
00661     /* The array will be packed as follows:
00662      * [numStrs, offset1, ... offsetN, characters data] 
00663      */
00664 
00665     bigArray.resize(totalLength);
00666 
00667     int* header = reinterpret_cast<int*>( &(bigArray[0]) );
00668     header[0] = x.size();
00669     for (Array<std::string>::size_type i=0; i<=x.size(); i++)
00670       {
00671         header[i+1] = offsets[i];
00672       }
00673 
00674     for (int i=0; i<x.length(); i++)
00675       {
00676         for (unsigned int j=0; j<x[i].length(); j++)
00677           {
00678             bigArray[offsets[i]+j] = x[i][j];
00679           }
00680       }
00681   }
00682 
00683   inline void MPIContainerComm<std::string>::unpack(const Array<char>& packed,
00684                                              Array<std::string>& x)
00685   {
00686     const int* header = reinterpret_cast<const int*>( &(packed[0]) );
00687 
00688     x.resize(header[0]);
00689     Array<int> offsets(x.size()+1);
00690     for (Array<std::string>::size_type i=0; i<=x.size(); i++) offsets[i] = header[i+1];
00691 
00692     for (Array<std::string>::size_type i=0; i<x.size(); i++)
00693       {
00694         x[i].resize(offsets[i+1]-offsets[i]);
00695         for (std::string::size_type j=0; j<x[i].length(); j++)
00696           {
00697             x[i][j] = packed[offsets[i] + j];
00698           }
00699       }
00700   }
00701 
00702   inline void MPIContainerComm<std::string>::getStrings(const Array<char>& bigArray,
00703                                                    const Array<int>& offsets,
00704                                                    Array<std::string>& x)
00705   {
00706     x.resize(offsets.length()-1);
00707     for (int i=0; i<x.length(); i++)
00708       {
00709         x[i].resize(offsets[i+1]-offsets[i]);
00710         for (unsigned int j=0; j<x[i].length(); j++)
00711           {
00712             x[i][j] = bigArray[offsets[i] + j];
00713           }
00714       }
00715   }
00716 #endif // DOXYGEN_SHOULD_SKIP_THIS
00717 
00718 }
00719 
00720 
00721 #endif
00722 
00723 

Site Contact