vtf-logo

pico/pico/DomainCoupler.h

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 //  ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00004 //
00005 //                                   Fehmi Cirak
00006 //                        California Institute of Technology
00007 //                           (C) 2003 All Rights Reserved
00008 //
00009 //  <LicenseText>
00010 //
00011 //  ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00012 //
00013 
00014 #ifndef _DomainCoupler_h_
00015 #define _DomainCoupler_h_
00016 
00017 #include <mpi.h>
00018 
00019 #include "PicoNames.h"
00020 #include "SubDomain.h"
00021 
00022 #include <map>
00023 #include <cassert>
00024 
00025 
00026 namespace pico {
00027     class DomainCoupler;
00028 }
00029 
00030 
00031 class pico::DomainCoupler{
00032     //
00033     // Each rank of the MPI communicator _pComm owns upto two types subdomains (0/1).
00034     //
00035     // Communication between subdomains of type 0 and 1 is called inter communication. 
00036     // Communication between subdomains of only  type 0 or only type 1 is called  
00037     // intra communication. 
00038     //
00039     
00040 public: 
00041     DomainCoupler(MPI_Comm pComm):_pComm(pComm){assert(_pComm != MPI_COMM_NULL);}
00042     ~DomainCoupler(){}
00043 
00044     // registration
00045     void registerSubDomain(const pico::SolvType& type, const double lo[3], 
00046                            const double up[3], const double& maxElemSize);
00047     void registerSubDomain(const pico::SolvType& type, const double * const coor, 
00048                            const int& nodes, const double& maxElemSize);
00049    
00050     // communication
00051     void prepareCommunicationsIntra(const pico::SolvType& type);
00052     void exchangeSubdomainDataIntra(const pico::SolvType& type);    
00053     void prepareCommunicationsInter();
00054     
00055     // register buffers
00056     template <typename T>
00057     void registerSendBuffer(const pico::SolvType& type, T* const start, int size);
00058 
00059     // accessors for receive buffer
00060     template <typename RT>
00061     RT* recvBuffer(const pico::SolvType& type, const RT& notUsed);
00062     template <typename T>
00063     int recvBufferSize(const pico::SolvType& type);
00064     
00065     // helper - for debugging purposes
00066     template <typename Iterator> 
00067     void relevantSubDomainNumbers(Iterator domainNumbers, int t, char *tch);
00068         
00069 // copy and equality constructor
00070 private:
00071     DomainCoupler(const DomainCoupler &);
00072     const DomainCoupler & operator=(const DomainCoupler &);
00073 
00074 private:
00075     void gatherSubDomains();
00076     void relevantSubDomainsInter();
00077     void relevantSubDomainsIntra(const pico::SolvType& type);
00078     void exchangeBufferSizesIntra(const pico::SolvType& type);
00079     
00080 private:
00081     typedef std::map<int, SubDomain *>       _SubDomainCont;
00082     typedef _SubDomainCont::iterator         _SubDomainIt;
00083   
00084     // std::map<processor number, subdomain pointer>
00085     _SubDomainCont                        _subDomains[2];
00086     
00087     _SubDomainCont                        _relevSubDomainsInter[2];
00088     _SubDomainCont                        _relevSubDomainsIntra[2];
00089      
00090     MPI_Comm                              _pComm;
00091 }; 
00092 
00093 
00094 
00095 namespace pico {
00096     // helper - for debugging purposes
00097     template <typename Iterator> 
00098     void DomainCoupler::relevantSubDomainNumbers(Iterator domainNumbers, int t, char *tch) 
00099     {
00100         std::string commT(tch);
00101         _SubDomainIt it, ite;
00102         
00103         if (commT== "Inter") {
00104             it = _relevSubDomainsInter[t].begin();
00105             ite = _relevSubDomainsInter[t].end();       
00106         } else if (commT=="Intra") {
00107             it = _relevSubDomainsIntra[t].begin();
00108             ite = _relevSubDomainsIntra[t].end();
00109         } else 
00110             return;
00111         
00112         for ( ; it!=ite; ++it) *(domainNumbers++) = it->first;
00113 
00114         return;
00115     }
00116             
00117 
00118     template <typename T>
00119     void DomainCoupler::registerSendBuffer(const pico::SolvType& type, 
00120                                            T* const start, int size) 
00121     {
00122         int rank;
00123         MPI_Comm_rank(_pComm, &rank);
00124         
00125         _SubDomainIt it, ite=_subDomains[type].end();
00126         
00127         it = _subDomains[type].find(rank);
00128         if (it==ite) return;
00129         
00130         SubDomain* myDomain = it->second;
00131         
00132         size_t unit = sizeof(T);
00133         size *= unit;
00134         
00135         char *startCh = reinterpret_cast<char *>(start);
00136 
00137         myDomain->registerSendBuffer(startCh, size);
00138         
00139         return;
00140     }
00141 
00142 
00143     template <typename RT>
00144     RT* DomainCoupler::recvBuffer(const pico::SolvType& type, const RT& notUsed) 
00145     {
00146         // notUsed is only "used" by the compiler
00147         char *tmp=NULL;
00148         int rank;
00149         MPI_Comm_rank(_pComm, &rank);
00150         
00151         _SubDomainIt it, ite=_subDomains[type].end();
00152         it = _subDomains[type].find(rank);
00153         if (it==ite) return NULL;
00154         SubDomain* myDomain = it->second;
00155         
00156         tmp = myDomain->recvBuffer();
00157         
00158         RT *tmpRT = reinterpret_cast<RT *>(tmp);
00159 
00160         return tmpRT;
00161     }
00162     
00163     
00164     template <typename T>
00165     int DomainCoupler::recvBufferSize(const pico::SolvType& type)  
00166     {
00167         int tmp=0, rank;
00168         MPI_Comm_rank(_pComm, &rank);
00169         
00170         _SubDomainIt it, ite=_subDomains[type].end();
00171         it = _subDomains[type].find(rank);
00172         if (it==ite) return tmp;
00173         SubDomain* myDomain = it->second;
00174         
00175         tmp = myDomain->recvBufferSize();
00176         
00177         size_t unit = sizeof(T);
00178         tmp /= unit;
00179 
00180         return tmp;
00181     }
00182 }
00183 
00184 
00185 #endif
00186 
00187 //  End of file

Generated on Fri Aug 24 13:00:24 2007 for SFC Thin-Shell Finite Element Solver by  doxygen 1.4.7