IFPACK  Development
 All Classes Files Functions Variables Enumerations Friends
Ifpack_OverlappingRowMatrix.cpp
00001 /*@HEADER
00002 // ***********************************************************************
00003 //
00004 //       Ifpack: Object-Oriented Algebraic Preconditioner Package
00005 //                 Copyright (2002) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
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 Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ***********************************************************************
00040 //@HEADER
00041 */
00042 
00043 #include "Ifpack_ConfigDefs.h"
00044 #include "Ifpack_OverlappingRowMatrix.h"
00045 #include "Epetra_RowMatrix.h"
00046 #include "Epetra_Map.h"
00047 #include "Epetra_CrsMatrix.h"
00048 #include "Epetra_Comm.h"
00049 #include "Epetra_MultiVector.h"
00050 
00051 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
00052 #include "Epetra_IntVector.h"
00053 #include "Epetra_MpiComm.h"
00054 #include "Teuchos_Hashtable.hpp"
00055 #include "Teuchos_Array.hpp"
00056 #include "EpetraExt_OperatorOut.h"
00057 #include "EpetraExt_BlockMapOut.h"
00058 #else
00059 # ifdef IFPACK_NODE_AWARE_CODE
00060 # include "Epetra_IntVector.h"
00061 # include "Epetra_MpiComm.h"
00062 # include "Teuchos_Hashtable.hpp"
00063 # include "Teuchos_Array.hpp"
00064 # include "EpetraExt_OperatorOut.h"
00065 # include "EpetraExt_BlockMapOut.h"
00066   extern int ML_NODE_ID;
00067   int IFPACK_NODE_ID;
00068 # endif
00069 #endif
00070 
00071 using namespace Teuchos;
00072 
00073 // ====================================================================== 
00074 
00075 template<typename int_type>
00076 void Ifpack_OverlappingRowMatrix::BuildMap(int OverlapLevel_in)
00077 {
00078   RCP<Epetra_Map> TmpMap;
00079   RCP<Epetra_CrsMatrix> TmpMatrix; 
00080   RCP<Epetra_Import> TmpImporter;
00081 
00082   // importing rows corresponding to elements that are 
00083   // in ColMap, but not in RowMap 
00084   const Epetra_Map *RowMap; 
00085   const Epetra_Map *ColMap; 
00086 
00087   std::vector<int_type> ExtElements; 
00088 
00089   for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap) {
00090     if (TmpMatrix != Teuchos::null) {
00091       RowMap = &(TmpMatrix->RowMatrixRowMap()); 
00092       ColMap = &(TmpMatrix->RowMatrixColMap()); 
00093     }
00094     else {
00095       RowMap = &(A().RowMatrixRowMap()); 
00096       ColMap = &(A().RowMatrixColMap()); 
00097     }
00098 
00099     int size = ColMap->NumMyElements() - RowMap->NumMyElements(); 
00100     std::vector<int_type> list(size); 
00101 
00102     int count = 0; 
00103 
00104     // define the set of rows that are in ColMap but not in RowMap 
00105     for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) { 
00106       int_type GID = (int_type) ColMap->GID64(i); 
00107       if (A().RowMatrixRowMap().LID(GID) == -1) { 
00108         typename std::vector<int_type>::iterator pos 
00109           = find(ExtElements.begin(),ExtElements.end(),GID); 
00110         if (pos == ExtElements.end()) { 
00111           ExtElements.push_back(GID);
00112           list[count] = GID; 
00113           ++count; 
00114         } 
00115       } 
00116     } 
00117 
00118     const int_type *listptr = NULL;
00119     if ( ! list.empty() ) listptr = &list[0];
00120     TmpMap = rcp( new Epetra_Map(-1,count, listptr,0,Comm()) ); 
00121 
00122     TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) ); 
00123 
00124     TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) ); 
00125 
00126     TmpMatrix->Import(A(),*TmpImporter,Insert); 
00127     TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap); 
00128 
00129   }
00130 
00131   // build the map containing all the nodes (original
00132   // matrix + extended matrix)
00133   std::vector<int_type> list(NumMyRowsA_ + ExtElements.size());
00134   for (int i = 0 ; i < NumMyRowsA_ ; ++i)
00135     list[i] = (int_type) A().RowMatrixRowMap().GID64(i);
00136   for (int i = 0 ; i < (int)ExtElements.size() ; ++i)
00137     list[i + NumMyRowsA_] = ExtElements[i];
00138 
00139   const int_type *listptr = NULL;
00140   if ( ! list.empty() ) listptr = &list[0];
00141   {
00142     Map_ = rcp( new Epetra_Map((int_type) -1, NumMyRowsA_ + ExtElements.size(),
00143                                listptr, 0, Comm()) );
00144   }
00145 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
00146   colMap_ = &*Map_;
00147 #else
00148 # ifdef IFPACK_NODE_AWARE_CODE
00149   colMap_ = &*Map_;
00150 # endif
00151 #endif
00152   // now build the map corresponding to all the external nodes
00153   // (with respect to A().RowMatrixRowMap().
00154   {
00155     const int_type * extelsptr = NULL;
00156     if ( ! ExtElements.empty() ) extelsptr = &ExtElements[0];
00157     ExtMap_ = rcp( new Epetra_Map((int_type) -1,ExtElements.size(),
00158                                   extelsptr,0,A().Comm()) );
00159   }
00160 }
00161 
00162 // ====================================================================== 
00163 // Constructor for the case of one core per subdomain
00164 Ifpack_OverlappingRowMatrix::
00165 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00166                             int OverlapLevel_in)  :
00167   Matrix_(Matrix_in),
00168   OverlapLevel_(OverlapLevel_in)
00169 {
00170   // should not be here if no overlap
00171   if (OverlapLevel_in == 0)
00172     IFPACK_CHK_ERRV(-1);
00173 
00174   // nothing to do as well with one process
00175   if (Comm().NumProc() == 1)
00176     IFPACK_CHK_ERRV(-1);
00177   
00178   NumMyRowsA_ = A().NumMyRows();
00179 
00180   // construct the external matrix
00181 
00182 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00183     if(A().RowMatrixRowMap().GlobalIndicesInt()) {
00184       BuildMap<int>(OverlapLevel_in);
00185     }
00186     else
00187 #endif
00188 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00189     if(A().RowMatrixRowMap().GlobalIndicesLongLong()) {
00190       BuildMap<long long>(OverlapLevel_in);
00191     }
00192     else
00193 #endif
00194       throw "Ifpack_OverlappingRowMatrix::Ifpack_OverlappingRowMatrix: GlobalIndices type unknown";
00195 
00196   ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*Map_,0) ); 
00197 
00198   ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) ); 
00199   ExtMatrix_->Import(A(),*ExtImporter_,Insert); 
00200   ExtMatrix_->FillComplete(A().OperatorDomainMap(),*Map_);
00201 
00202   Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) );
00203 
00204   // fix indices for overlapping matrix
00205   NumMyRowsB_ = B().NumMyRows();
00206   NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;
00207   NumMyCols_ = NumMyRows_;
00208   
00209   NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
00210   
00211   NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
00212   long long NumMyNonzeros_tmp = NumMyNonzeros_;
00213   Comm().SumAll(&NumMyNonzeros_tmp,&NumGlobalNonzeros_,1);
00214   MaxNumEntries_ = A().MaxNumEntries();
00215   
00216   if (MaxNumEntries_ < B().MaxNumEntries())
00217     MaxNumEntries_ = B().MaxNumEntries();
00218 
00219 }
00220 
00221 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
00222 
00223 // ====================================================================== 
00224 // Constructor for the case of two or more cores per subdomain
00225 Ifpack_OverlappingRowMatrix::
00226 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00227                             int OverlapLevel_in, int subdomainID)  :
00228   Matrix_(Matrix_in),
00229   OverlapLevel_(OverlapLevel_in)
00230 {
00231 
00232   //FIXME timing
00233 #ifdef IFPACK_OVA_TIME_BUILD
00234   double t0 = MPI_Wtime();
00235   double t1, true_t0=t0;
00236 #endif
00237   //FIXME timing
00238 
00239   const int votesMax = INT_MAX;
00240 
00241   // should not be here if no overlap
00242   if (OverlapLevel_in == 0)
00243     IFPACK_CHK_ERRV(-1);
00244 
00245   // nothing to do as well with one process
00246   if (Comm().NumProc() == 1)
00247     IFPACK_CHK_ERRV(-1);
00248 
00249   // subdomainID is the node (aka socket) number, and so is system dependent
00250   // nodeComm is the communicator for all the processes on a particular node
00251   // these processes will have the same subdomainID.
00252 # ifdef HAVE_MPI
00253   const Epetra_MpiComm *pComm = dynamic_cast<const Epetra_MpiComm*>( &Comm() );
00254   assert(pComm != NULL);
00255   MPI_Comm nodeMPIComm;
00256   MPI_Comm_split(pComm->Comm(),subdomainID,pComm->MyPID(),&nodeMPIComm);
00257   Epetra_MpiComm *nodeComm = new Epetra_MpiComm(nodeMPIComm);
00258 # else
00259   Epetra_SerialComm *nodeComm =  dynamic_cast<Epetra_MpiComm*>( &(Matrix_in->RowMatrixRowMap().Comm()) );
00260 # endif
00261   
00262   NumMyRowsA_ = A().NumMyRows();
00263 
00264   RCP<Epetra_Map> TmpMap;
00265   RCP<Epetra_CrsMatrix> TmpMatrix; 
00266   RCP<Epetra_Import> TmpImporter;
00267 
00268   // importing rows corresponding to elements that are 
00269   // in ColMap, but not in RowMap 
00270   const Epetra_Map *RowMap; 
00271   const Epetra_Map *ColMap; 
00272   const Epetra_Map *DomainMap;
00273 
00274   // TODO Count #connections from nodes I own to each ghost node
00275 
00276 
00277   //FIXME timing
00278 #ifdef IFPACK_OVA_TIME_BUILD
00279   t1 = MPI_Wtime(); 
00280   fprintf(stderr,"[node %d]: time for initialization %2.3e\n", subdomainID, t1-t0);
00281   t0=t1;
00282 #endif
00283   //FIXME timing
00284 
00285 #ifdef IFPACK_OVA_TIME_BUILD
00286   t1 = MPI_Wtime();
00287   fprintf(stderr,"[node %d]: overlap hash table allocs %2.3e\n", subdomainID, t1-t0);
00288   t0=t1;
00289 #endif
00290 
00291   Teuchos::Hashtable<int,int> colMapTable(3 * A().RowMatrixColMap().NumMyElements() );
00292 
00293   // ghostTable holds off-node GIDs that are connected to on-node rows and can potentially be this PID's overlap
00294   // TODO hopefully 3 times the # column entries is a reasonable table size
00295   Teuchos::Hashtable<int,int> ghostTable(3 * A().RowMatrixColMap().NumMyElements() );
00296 
00297   /* ** ************************************************************************** ** */
00298   /* ** ********************** start of main overlap loop ************************ ** */
00299   /* ** ************************************************************************** ** */
00300   for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap)
00301   {
00302     // nbTable holds node buddy GIDs that are connected to current PID's rows, i.e., GIDs that should be in the overlapped
00303     // matrix's column map
00304 
00305     if (TmpMatrix != Teuchos::null) {
00306       RowMap = &(TmpMatrix->RowMatrixRowMap()); 
00307       ColMap = &(TmpMatrix->RowMatrixColMap()); 
00308       DomainMap = &(TmpMatrix->OperatorDomainMap());
00309     }
00310     else {
00311       RowMap = &(A().RowMatrixRowMap()); 
00312       ColMap = &(A().RowMatrixColMap()); 
00313       DomainMap = &(A().OperatorDomainMap());
00314     }
00315 
00316 #ifdef IFPACK_OVA_TIME_BUILD
00317     t1 = MPI_Wtime();
00318     fprintf(stderr,"[node %d]: overlap pointer pulls %2.3e\n", subdomainID, t1-t0);
00319     t0=t1;
00320 #endif
00321     
00322     // For each column ID, determine the owning node (as opposed to core)
00323     // ID of the corresponding row.
00324     Epetra_IntVector colIdList( *ColMap );
00325     Epetra_IntVector rowIdList(*DomainMap);
00326     rowIdList.PutValue(subdomainID);  
00327     Teuchos::RCP<Epetra_Import> nodeIdImporter = rcp(new Epetra_Import( *ColMap, *DomainMap ));
00328 
00329 #ifdef IFPACK_OVA_TIME_BUILD
00330     t1 = MPI_Wtime();
00331     fprintf(stderr,"[node %d]: overlap intvector/importer allocs %2.3e\n", subdomainID, t1-t0);
00332     t0=t1;
00333 #endif
00334     
00335     colIdList.Import(rowIdList,*nodeIdImporter,Insert);
00336     
00337     int size = ColMap->NumMyElements() - RowMap->NumMyElements(); 
00338     int count = 0; 
00339 
00340 #ifdef IFPACK_OVA_TIME_BUILD
00341     t1 = MPI_Wtime();
00342     fprintf(stderr,"[node %d]: overlap col/row ID imports %2.3e\n", subdomainID, t1-t0);
00343     t0=t1;
00344 #endif
00345 
00346     
00347     // define the set of off-node rows that are in ColMap but not in RowMap
00348     // This naturally does not include off-core rows that are on the same node as me, i.e., node buddy rows.
00349     for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) {
00350       int GID = ColMap->GID(i); 
00351       int myvotes=0;
00352       if (ghostTable.containsKey(GID)) myvotes = ghostTable.get(GID);
00353       if ( colIdList[i] != subdomainID && myvotes < votesMax) // don't include anybody found in a previous round
00354       {
00355         int votes;
00356         if (ghostTable.containsKey(GID)) {
00357           votes = ghostTable.get(GID);
00358           votes++;
00359           ghostTable.put(GID,votes);
00360         } else {
00361           ghostTable.put(GID,1);
00362         }
00363       }
00364     } //for (int i = 0 ; i < ColMap->NumMyElements() ; ++i)
00365 
00366     Teuchos::Array<int> gidsarray,votesarray;
00367     ghostTable.arrayify(gidsarray,votesarray);
00368     int *gids = gidsarray.getRawPtr();
00369     int *votes = votesarray.getRawPtr();
00370 
00371     /*
00372        This next bit of code decides which node buddy (NB) gets which ghost points.  Everyone sends their
00373        list of ghost points to pid 0 of the local subcommunicator.  Pid 0 decides who gets what:
00374 
00375           - if a ghost point is touched by only one NB, that NB gets the ghost point
00376           - if two or more NBs share a ghost point, whichever NB has the most connections to the ghost
00377             point gets it.
00378     */
00379 
00380 #   ifdef HAVE_MPI  //FIXME What if we build in serial?!  This file will likely not compile.
00381     int *cullList;
00382     int ncull;
00383     //int mypid = nodeComm->MyPID();
00384 
00385     //FIXME timing
00386 #ifdef IFPACK_OVA_TIME_BUILD
00387     t1 = MPI_Wtime();
00388     fprintf(stderr,"[node %d]: overlap before culling %2.3e\n", subdomainID, t1-t0);
00389     t0=t1;
00390 #endif
00391     //FIXME timing
00392 
00393 
00394     if (nodeComm->MyPID() == 0)
00395     {
00396       // Figure out how much pid 0 is to receive
00397       MPI_Status status;
00398       int *lengths = new int[nodeComm->NumProc()+1];
00399       lengths[0] = 0;
00400       lengths[1] = ghostTable.size();
00401       for (int i=1; i<nodeComm->NumProc(); i++) {
00402         int leng;
00403         MPI_Recv( &leng, 1, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00404         lengths[i+1] = lengths[i] + leng;
00405       }
00406       int total = lengths[nodeComm->NumProc()];
00407 
00408       int* ghosts = new int[total];
00409       for (int i=0; i<total; i++) ghosts[i] = -9;
00410       int *round  = new int[total];
00411       int *owningpid  = new int[total];
00412 
00413       for (int i=1; i<nodeComm->NumProc(); i++) {
00414         int count = lengths[i+1] - lengths[i];
00415         MPI_Recv( ghosts+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00416         MPI_Recv( round+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00417       }
00418 
00419       // put in pid 0's info
00420       for (int i=0; i<lengths[1]; i++) {
00421         ghosts[i] = gids[i];
00422         round[i] = votes[i];
00423         owningpid[i] = 0;
00424       }
00425 
00426       // put in the pid associated with each ghost
00427       for (int j=1; j<nodeComm->NumProc(); j++) {
00428         for (int i=lengths[j]; i<lengths[j+1]; i++) {
00429           owningpid[i] = j;
00430         }
00431       }
00432 
00433       // sort everything based on the ghost gids
00434       int* roundpid[2];
00435       roundpid[0] = round; roundpid[1] = owningpid;
00436       Epetra_Util epetraUtil;
00437       epetraUtil.Sort(true,total,ghosts,0,0,2,roundpid);
00438 
00439       //set up arrays that get sent back to node buddies and that tell them what ghosts to cull
00440       int *nlosers = new int[nodeComm->NumProc()];
00441       int** losers = new int*[nodeComm->NumProc()];
00442       for (int i=0; i<nodeComm->NumProc(); i++) {
00443         nlosers[i] = 0;
00444         losers[i] = new int[ lengths[i+1]-lengths[i] ];
00445       }
00446 
00447       // Walk through ghosts array and and for each sequence of duplicate ghost GIDs, choose just one NB to keep it.
00448       // The logic is pretty simple.  The ghost list is sorted, so all duplicate PIDs are together.
00449       // The list is traversed.  As duplicates are found, node pid 0 keeps track of the current "winning"
00450       // pid.  When a pid is determined to have "lost" (less votes/connections to the current GID), the
00451       // GID is added to that pid's list of GIDs to be culled.  At the end of the repeated sequence, we have
00452       // a winner, and other NBs know whether they need to delete it from their import list.
00453       int max=0;   //for duplicated ghosts, index of pid with most votes and who hence keeps the ghost.
00454                    // TODO to break ties randomly
00455 
00456       for (int i=1; i<total; i++) {
00457         if (ghosts[i] == ghosts[i-1]) {
00458           int idx = i; // pid associated with idx is current "loser"
00459           if (round[i] > round[max]) {
00460             idx = max;
00461             max=i;
00462           }
00463           int j = owningpid[idx];
00464           losers[j][nlosers[j]++] = ghosts[idx];
00465         } else {
00466           max=i;
00467         }
00468       } //for (int i=1; i<total; i++)
00469 
00470       delete [] round;
00471       delete [] ghosts;
00472       delete [] owningpid;
00473 
00474       // send the arrays of ghost GIDs to be culled back to the respective node buddies
00475       for (int i=1; i<nodeComm->NumProc(); i++) {
00476         MPI_Send( nlosers+i, 1, MPI_INT, i, 8675, nodeComm->Comm());
00477         MPI_Send( losers[i], nlosers[i], MPI_INT, i, 8675, nodeComm->Comm());
00478       }
00479 
00480       //FIXME Unnecessary memory allocation and copying, but makes culling code cleaner
00481       //Could we stick this info into gids and votes, since neither is being used anymore?
00482       //TODO Instead of using "losers" arrays, just use "cullList" as in the else clause
00483       ncull = nlosers[0];
00484       cullList = new int[ncull+1];
00485       for (int i=0; i<nlosers[0]; i++)
00486         cullList[i] = losers[0][i];
00487 
00488       for (int i=0; i<nodeComm->NumProc(); i++)
00489         delete [] losers[i];
00490 
00491       delete [] lengths;
00492       delete [] losers;
00493       delete [] nlosers;
00494 
00495     } else { //everyone but pid 0
00496 
00497       // send to node pid 0 all ghosts that this pid could potentially import
00498       int hashsize = ghostTable.size();
00499       MPI_Send( &hashsize, 1, MPI_INT, 0, 8675, nodeComm->Comm());
00500       MPI_Send( gids, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
00501       MPI_Send( votes, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
00502 
00503       // receive the ghost GIDs that should not be imported (subset of the list sent off just above)
00504       MPI_Status status;
00505       MPI_Recv( &ncull, 1, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
00506       cullList = new int[ncull+1];
00507       MPI_Recv( cullList, ncull, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
00508     }
00509 
00510 
00511     //TODO clean out hash table after each time through overlap loop   4/1/07 JJH done moved both hash tables to local scope
00512 
00513     // Remove from my row map all off-node ghosts that will be imported by a node buddy.
00514     for (int i=0; i<ncull; i++) {
00515       try{ghostTable.remove(cullList[i]);}
00516 
00517       catch(...) {
00518         printf("pid %d: In OverlappingRowMatrix ctr, problem removing ghost elt %d from ghostTable\n",
00519                Comm().MyPID(),cullList[i]);
00520         fflush(stdout);
00521       }
00522     }//for
00523 
00524     delete [] cullList;
00525 
00526     // Save off the remaining ghost GIDs from the current overlap round.
00527     // These are off-node GIDs (rows) that I will import.
00528     gidsarray.clear(); votesarray.clear();
00529     ghostTable.arrayify(gidsarray,votesarray);
00530 
00531     std::vector<int> list(size); 
00532     count=0;
00533     for (int i=0; i<ghostTable.size(); i++) {
00534       // if votesarray[i] == votesmax, then this GID was found during a previous overlap round
00535       if (votesarray[i] < votesMax) {
00536         list[count] = gidsarray[i];
00537         ghostTable.put(gidsarray[i],votesMax); //set this GID's vote to the maximum so that this PID alway wins in any subsequent overlap tiebreaking.
00538         count++;
00539       }
00540     }
00541 
00542     //FIXME timing
00543 #ifdef IFPACK_OVA_TIME_BUILD
00544     t1 = MPI_Wtime();
00545     fprintf(stderr,"[node %d]: overlap duplicate removal %2.3e\n", subdomainID, t1-t0);
00546     t0=t1;
00547 #endif
00548     //FIXME timing
00549 
00550 #   endif //ifdef HAVE_MPI
00551 
00552     TmpMap = rcp( new Epetra_Map(-1,count, &list[0],0,Comm()) );
00553 
00554     TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) ); 
00555 
00556     TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) ); 
00557 
00558     TmpMatrix->Import(A(),*TmpImporter,Insert); 
00559     TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap); 
00560 
00561     //FIXME timing
00562 #ifdef IFPACK_OVA_TIME_BUILD
00563     t1 = MPI_Wtime();
00564     fprintf(stderr,"[node %d]: overlap TmpMatrix fillcomplete %2.3e\n", subdomainID, t1-t0);
00565     t0=t1;
00566 #endif
00567     //FIXME timing
00568 
00569 
00570     // These next two imports get the GIDs that need to go into the column map of the overlapped matrix.
00571 
00572     // For each column ID in the overlap, determine the owning node (as opposed to core)
00573     // ID of the corresponding row.  Save those column IDs whose owning node is the current one.
00574     // This should get all the imported ghost GIDs.
00575     Epetra_IntVector ov_colIdList( TmpMatrix->ColMap() );
00576     ov_colIdList.PutValue(-1);
00577     Epetra_IntVector ov_rowIdList( TmpMatrix->RowMap() );
00578     ov_rowIdList.PutValue(subdomainID);  
00579     Teuchos::RCP<Epetra_Import> ov_nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), TmpMatrix->RowMap()));
00580     ov_colIdList.Import(ov_rowIdList,*ov_nodeIdImporter,Insert);
00581 
00582     for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
00583        if (ov_colIdList[i] == subdomainID) {
00584          int GID = (ov_colIdList.Map()).GID(i);
00585          colMapTable.put(GID,1);
00586        }
00587     }
00588 
00589     // Do a second import of the owning node ID from A's rowmap to TmpMat's column map.  This ensures that
00590     // all GIDs that belong to a node buddy and are in a ghost row's sparsity pattern will be in the final
00591     // overlapped matrix's column map.
00592     ov_colIdList.PutValue(-1);
00593     Epetra_IntVector ArowIdList( A().RowMatrixRowMap() );
00594     ArowIdList.PutValue(subdomainID);
00595     nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), A().RowMatrixRowMap() ));
00596     ov_colIdList.Import(ArowIdList,*nodeIdImporter,Insert);
00597 
00598     for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
00599        if (ov_colIdList[i] == subdomainID) {
00600          int GID = (ov_colIdList.Map()).GID(i);
00601          colMapTable.put(GID,1);
00602        }
00603     }
00604 
00605     //FIXME timing
00606 #ifdef IFPACK_OVA_TIME_BUILD
00607     t1 = MPI_Wtime();
00608     fprintf(stderr,"[node %d]: overlap 2 imports to fix up colmap %2.3e\n", subdomainID, t1-t0);
00609     t0=t1;
00610 #endif
00611     //FIXME timing
00612 
00613   } //for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap)
00614 
00615   /* ** ************************************************************************ ** */
00616   /* ** ********************** end of main overlap loop ************************ ** */
00617   /* ** ************************************************************************ ** */
00618 
00619   // off-node GIDs that will be in the overlap
00620   std::vector<int> ghostElements; 
00621 
00622   Teuchos::Array<int> gidsarray,votesarray;
00623   ghostTable.arrayify(gidsarray,votesarray);
00624   for (int i=0; i<ghostTable.size(); i++) {
00625     ghostElements.push_back(gidsarray[i]);
00626   }
00627 
00628     for (int i = 0 ; i < A().RowMatrixColMap().NumMyElements() ; ++i) {
00629       int GID = A().RowMatrixColMap().GID(i);
00630       // Remove any entries that are in A's original column map
00631       if (colMapTable.containsKey(GID)) {
00632         try{colMapTable.remove(GID);}
00633         catch(...) {
00634           printf("pid %d: In OverlappingRowMatrix ctr, problem removing entry %d from colMapTable\n", Comm().MyPID(),GID);
00635           fflush(stdout);
00636         }
00637       }
00638     }
00639 
00640     // GIDs that will be in the overlapped matrix's column map
00641     std::vector<int> colMapElements; 
00642 
00643   gidsarray.clear(); votesarray.clear();
00644   colMapTable.arrayify(gidsarray,votesarray);
00645   for (int i=0; i<colMapTable.size(); i++)
00646     colMapElements.push_back(gidsarray[i]);
00647 
00648 /*
00649    We need two maps here.  The first is the row map, which we've got by using my original rows
00650    plus whatever I've picked up in ghostElements.
00651 
00652    The second is a column map.  This map should include my rows, plus any columns that could come from node buddies.
00653    These GIDs come from the std:array colMapElements, which in turn comes from colMapTable.
00654    This map should *not* omit ghosts that have been imported by my node buddies, i.e., for any row that I own,
00655    the stencil should include all column GIDs (including imported ghosts) that are on the node.
00656 */
00657 
00658   // build the row map containing all the nodes (original matrix + off-node matrix)
00659   std::vector<int> rowList(NumMyRowsA_ + ghostElements.size());
00660   for (int i = 0 ; i < NumMyRowsA_ ; ++i)
00661     rowList[i] = A().RowMatrixRowMap().GID(i);
00662   for (int i = 0 ; i < (int)ghostElements.size() ; ++i)
00663     rowList[i + NumMyRowsA_] = ghostElements[i];
00664 
00665   // row map for the overlapped matrix (local + overlap)
00666   Map_ = rcp( new Epetra_Map(-1, NumMyRowsA_ + ghostElements.size(), &rowList[0], 0, Comm()) );
00667 
00668   // build the column map for the overlapping matrix
00669   //vector<int> colList(colMapElements.size());
00670   // column map for the overlapped matrix (local + overlap)
00671   //colMap_ = rcp( new Epetra_Map(-1, colMapElements.size(), &colList[0], 0, Comm()) );
00672   //for (int i = 0 ; i < (int)colMapElements.size() ; i++)
00673   //  colList[i] = colMapElements[i];
00674   std::vector<int> colList(A().RowMatrixColMap().NumMyElements() + colMapElements.size());
00675   int nc = A().RowMatrixColMap().NumMyElements();
00676   for (int i = 0 ; i < nc; i++)
00677     colList[i] = A().RowMatrixColMap().GID(i);
00678   for (int i = 0 ; i < (int)colMapElements.size() ; i++)
00679     colList[nc+i] = colMapElements[i];
00680 
00681   // column map for the overlapped matrix (local + overlap)
00682   //colMap_ = rcp( new Epetra_Map(-1, A().RowMatrixColMap().NumMyElements() + colMapElements.size(), &colList[0], 0, Comm()) );
00683   colMap_ = new Epetra_Map(-1, A().RowMatrixColMap().NumMyElements() + colMapElements.size(), &colList[0], 0, Comm());
00684     
00685 
00686   //FIXME timing
00687 #ifdef IFPACK_OVA_TIME_BUILD
00688   t1 = MPI_Wtime();
00689   fprintf(stderr,"[node %d]: time to init B() col/row maps %2.3e\n", subdomainID, t1-t0);
00690   t0=t1;
00691 #endif
00692   //FIXME timing
00693 
00694   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00695   //++++ start of sort
00696   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00697   // build the column map, but don't use a copy constructor, b/c local communicator SubComm_ is
00698   // different from that of Matrix.
00699   try {
00700     // build row map based on the local communicator.  We need this temporarily to build the column map.
00701     Epetra_Map* nodeMap = new Epetra_Map(-1,NumMyRowsA_ + ghostElements.size(),&rowList[0],0,*nodeComm);
00702     int numMyElts = colMap_->NumMyElements();
00703     assert(numMyElts!=0);
00704 
00705     // The column map *must* be sorted: first locals, then ghosts.  The ghosts
00706     // must be further sorted so that they are contiguous by owning processor.
00707 
00708     // For each GID in column map, determine owning PID in local communicator.
00709     int* myGlobalElts = new int[numMyElts];
00710     colMap_->MyGlobalElements(myGlobalElts);
00711     int* pidList = new int[numMyElts];
00712     nodeMap->RemoteIDList(numMyElts, myGlobalElts, pidList, 0);
00713 
00714     // First sort column map GIDs based on the owning PID in local communicator.
00715     Epetra_Util Util;
00716     int *tt[1];
00717     tt[0] = myGlobalElts;
00718     Util.Sort(true, numMyElts, pidList, 0, (double**)0, 1, tt);
00719 
00720     // For each remote PID, sort the column map GIDs in ascending order.
00721     // Don't sort the local PID's entries.
00722     int localStart=0;
00723     while (localStart<numMyElts) {
00724       int currPID = (pidList)[localStart];
00725       int i=localStart;
00726       while (i<numMyElts && (pidList)[i] == currPID) i++;
00727       if (currPID != nodeComm->MyPID())
00728         Util.Sort(true, i-localStart, (myGlobalElts)+localStart, 0, 0, 0, 0);
00729       localStart = i;
00730     }
00731 
00732     // now move the local entries to the front of the list
00733     localStart=0;
00734     while (localStart<numMyElts && (pidList)[localStart] != nodeComm->MyPID()) localStart++;
00735     assert(localStart != numMyElts);
00736     int localEnd=localStart;
00737     while (localEnd<numMyElts && (pidList)[localEnd] == nodeComm->MyPID()) localEnd++;
00738     int* mySortedGlobalElts = new int[numMyElts];
00739     int leng = localEnd - localStart;
00740     /* This is a little gotcha.  It appears that the ordering of the column map's local entries
00741        must be the same as that of the domain map's local entries.  See the comment in method
00742        MakeColMap() in Epetra_CrsGraph.cpp, line 1072. */
00743     int *rowGlobalElts =  nodeMap->MyGlobalElements();
00744 
00745     /*TODO TODO TODO NTS rows 68 and 83 show up as local GIDs in rowGlobalElts for both pids 0 & 1 on node 0.
00746                     This seems to imply that there is something wrong with rowList!!
00747                     It is ok for the overlapped matrix's row map to have repeated entries (it's overlapped, after all).
00748                     But ... on a node, there should be no repeats.
00749                     Ok, here's the underlying cause.  On node 0, gpid 1 finds 83 on overlap round 0.  gpid 0 finds 83
00750                     on overlap round 1.  This shouldn't happen .. once a node buddy "owns" a row, no other node buddy
00751                     should ever import ("own") that row.
00752 
00753                     Here's a possible fix.  Extend tie breaking across overlap rounds.  This means sending to lpid 0
00754                     a list of *all* rows you've imported (this round plus previous rounds) for tie-breaking
00755                     If a nb won a ghost row in a previous round, his votes for that ghost row should be really high, i.e.,
00756                     he should probably win that ghost row forever.
00757                     7/17/09
00758 
00759       7/28/09 JJH I believe I've fixed this, but I thought it might be helpful to have this comment in here, in case I missed something.
00760     */
00761 
00762     //move locals to front of list
00763     for (int i=0; i<leng; i++) {
00764       /* //original code */
00765       (mySortedGlobalElts)[i] = rowGlobalElts[i];
00766       //(&*mySortedGlobalElts)[i] = (&*myGlobalElts)[localStart+i];
00767       //(&*mySortedPidList)[i] = (&*pidList)[localStart+i];
00768     }
00769     for (int i=leng; i< localEnd; i++) {
00770       (mySortedGlobalElts)[i] = (myGlobalElts)[i-leng];
00771     }
00772     for (int i=localEnd; i<numMyElts; i++) {
00773       (mySortedGlobalElts)[i] = (myGlobalElts)[i];
00774     }
00775 
00776     //FIXME timing
00777 #ifdef IFPACK_OVA_TIME_BUILD
00778     t1 = MPI_Wtime();
00779     fprintf(stderr,"[node %d]: time to sort col map arrays (cp=1) %2.3e\n", subdomainID, t1-t0);
00780     t0=t1;
00781 #endif
00782     //FIXME timing
00783 
00784     int indexBase = colMap_->IndexBase();
00785     if (colMap_) delete colMap_;
00786     colMap_ = new Epetra_Map(-1,numMyElts,mySortedGlobalElts,indexBase,Comm());
00787 
00788     delete nodeMap;
00789     delete [] myGlobalElts;
00790     delete [] pidList;
00791     delete [] mySortedGlobalElts;
00792 
00793   } //try
00794   catch(...) {
00795     printf("** * Ifpack_OverlappingRowmatrix ctor: problem creating column map * **\n\n");
00796   }
00797 
00798   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00799   //++++ end of sort
00800   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00801 
00802   //FIXME timing
00803 #ifdef IFPACK_OVA_TIME_BUILD
00804   t1 = MPI_Wtime();
00805   fprintf(stderr,"[node %d]: time to sort col map arrays (cp=2) %2.3e\n", subdomainID, t1-t0);
00806   t0=t1;
00807 #endif
00808   //FIXME timing
00809 
00810 
00811 /*
00812 
00813    FIXME
00814    Does the column map need to be sorted for the overlapping matrix?
00815 
00816    The column map *must* be sorted:
00817 
00818         first locals
00819         then ghosts
00820 
00821    The ghosts must be further sorted so that they are contiguous by owning processor
00822 
00823   int* RemoteSizeList = 0
00824   int* RemoteColIndices = ColIndices.Values() + NumLocalColGIDs; // Points to back end of ColIndices
00825 
00826   EPETRA_CHK_ERR(DomainMap.RemoteIDList(NumRemoteColGIDs, RemoteColIndices, PIDList.Values(), 0, RemoteSizeList));
00827   Epetra_Util epetraUtil;
00828   SortLists[0] = RemoteColIndices;
00829   SortLists[1] = RemoteSizeList;
00830   epetraUtil.Sort(true, NumRemoteColGIDs, PIDList.Values(), 0, 0, NLists, SortLists);
00831 */
00832 
00833   // now build the map corresponding to all the external nodes
00834   // (with respect to A().RowMatrixRowMap().
00835   ExtMap_ = rcp( new Epetra_Map(-1,ghostElements.size(), &ghostElements[0],0,Comm()) );
00836   ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*colMap_,0) ); 
00837 
00838   ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) ); 
00839   ExtMatrix_->Import(A(),*ExtImporter_,Insert); 
00840 
00841   //FIXME timing
00842 #ifdef IFPACK_OVA_TIME_BUILD
00843   t1 = MPI_Wtime();
00844   fprintf(stderr,"[node %d]: time to import and setup ExtMap_ %2.3e\n", subdomainID, t1-t0);
00845   t0=t1;
00846 #endif
00847   //FIXME timing
00848 
00849 
00850 /*
00851   Notes to self:    In FillComplete, the range map does not have to be 1-1 as long as
00852                     (row map == range map).  Ditto for the domain map not being 1-1
00853                     if (col map == domain map).
00854                     
00855 */
00856 
00857   ExtMatrix_->FillComplete( *colMap_ , *Map_ ); //FIXME wrong
00858 
00859   //FIXME timing
00860 #ifdef IFPACK_OVA_TIME_BUILD
00861   t1 = MPI_Wtime();
00862   fprintf(stderr,"[node %d]: time to FillComplete B() %2.3e\n", subdomainID, t1-t0);
00863   t0=t1;
00864 #endif
00865   //FIXME timing
00866 
00867 
00868   // Note: B() = *ExtMatrix_ .
00869 
00870   Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) ); //FIXME is this right?!
00871 
00872   // fix indices for overlapping matrix
00873   NumMyRowsB_ = B().NumMyRows();
00874   NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;  //TODO is this wrong for a subdomain on >1 processor? // should be ok
00875   //NumMyCols_ = NumMyRows_;  //TODO is this wrong for a subdomain on >1 processor?  // YES!!!
00876   //NumMyCols_ = A().NumMyCols() + B().NumMyCols();
00877   NumMyCols_ = B().NumMyCols();
00878 
00879   /*FIXME*/ //somehow B's NumMyCols is the entire subdomain (local + overlap)
00880 
00881   NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
00882   
00883   NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
00884   long long NumMyNonZerosTemp = NumMyNonzeros_;
00885   Comm().SumAll(&NumMyNonZerosTemp,&NumGlobalNonzeros_,1);
00886   MaxNumEntries_ = A().MaxNumEntries();
00887   
00888   if (MaxNumEntries_ < B().MaxNumEntries())
00889     MaxNumEntries_ = B().MaxNumEntries();
00890 
00891 # ifdef HAVE_MPI
00892   delete nodeComm;
00893 # endif
00894 
00895   //FIXME timing
00896 #ifdef IFPACK_OVA_TIME_BUILD
00897   t1 = MPI_Wtime();
00898   fprintf(stderr,"[node %d]: time for final calcs %2.3e\n", subdomainID, t1-t0);
00899   fprintf(stderr,"[node %d]: total IORM ctor time %2.3e\n", subdomainID, t1-true_t0);
00900 #endif
00901   //FIXME timing
00902 
00903 
00904 } //Ifpack_OverlappingRowMatrix() ctor for more than one core
00905 
00906 #else
00907 # ifdef IFPACK_NODE_AWARE_CODE
00908 
00909 // ====================================================================== 
00910 // Constructor for the case of two or more cores per subdomain
00911 Ifpack_OverlappingRowMatrix::
00912 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00913                             int OverlapLevel_in, int nodeID)  :
00914   Matrix_(Matrix_in),
00915   OverlapLevel_(OverlapLevel_in)
00916 {
00917 
00918   //FIXME timing
00919 #ifdef IFPACK_OVA_TIME_BUILD
00920   double t0 = MPI_Wtime();
00921   double t1, true_t0=t0;
00922 #endif
00923   //FIXME timing
00924 
00925   const int votesMax = INT_MAX;
00926 
00927   // should not be here if no overlap
00928   if (OverlapLevel_in == 0)
00929     IFPACK_CHK_ERRV(-1);
00930 
00931   // nothing to do as well with one process
00932   if (Comm().NumProc() == 1)
00933     IFPACK_CHK_ERRV(-1);
00934 
00935   // nodeID is the node (aka socket) number, and so is system dependent
00936   // nodeComm is the communicator for all the processes on a particular node
00937   // these processes will have the same nodeID.
00938 # ifdef HAVE_MPI
00939   const Epetra_MpiComm *pComm = dynamic_cast<const Epetra_MpiComm*>( &Comm() );
00940   assert(pComm != NULL);
00941   MPI_Comm nodeMPIComm;
00942   MPI_Comm_split(pComm->Comm(),nodeID,pComm->MyPID(),&nodeMPIComm);
00943   Epetra_MpiComm *nodeComm = new Epetra_MpiComm(nodeMPIComm);
00944 # else
00945   Epetra_SerialComm *nodeComm =  dynamic_cast<Epetra_MpiComm*>( &(Matrix_in->RowMatrixRowMap().Comm()) );
00946 # endif
00947   
00948   NumMyRowsA_ = A().NumMyRows();
00949 
00950   RCP<Epetra_Map> TmpMap;
00951   RCP<Epetra_CrsMatrix> TmpMatrix; 
00952   RCP<Epetra_Import> TmpImporter;
00953 
00954   // importing rows corresponding to elements that are 
00955   // in ColMap, but not in RowMap 
00956   const Epetra_Map *RowMap; 
00957   const Epetra_Map *ColMap; 
00958   const Epetra_Map *DomainMap;
00959 
00960   // TODO Count #connections from nodes I own to each ghost node
00961 
00962 
00963   //FIXME timing
00964 #ifdef IFPACK_OVA_TIME_BUILD
00965   t1 = MPI_Wtime(); 
00966   fprintf(stderr,"[node %d]: time for initialization %2.3e\n", nodeID, t1-t0);
00967   t0=t1;
00968 #endif
00969   //FIXME timing
00970 
00971 #ifdef IFPACK_OVA_TIME_BUILD
00972   t1 = MPI_Wtime();
00973   fprintf(stderr,"[node %d]: overlap hash table allocs %2.3e\n", nodeID, t1-t0);
00974   t0=t1;
00975 #endif
00976 
00977   Teuchos::Hashtable<int,int> colMapTable(3 * A().RowMatrixColMap().NumMyElements() );
00978 
00979   // ghostTable holds off-node GIDs that are connected to on-node rows and can potentially be this PID's overlap
00980   // TODO hopefully 3 times the # column entries is a reasonable table size
00981   Teuchos::Hashtable<int,int> ghostTable(3 * A().RowMatrixColMap().NumMyElements() );
00982 
00983   /* ** ************************************************************************** ** */
00984   /* ** ********************** start of main overlap loop ************************ ** */
00985   /* ** ************************************************************************** ** */
00986   for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap)
00987   {
00988     // nbTable holds node buddy GIDs that are connected to current PID's rows, i.e., GIDs that should be in the overlapped
00989     // matrix's column map
00990 
00991     if (TmpMatrix != Teuchos::null) {
00992       RowMap = &(TmpMatrix->RowMatrixRowMap()); 
00993       ColMap = &(TmpMatrix->RowMatrixColMap()); 
00994       DomainMap = &(TmpMatrix->OperatorDomainMap());
00995     }
00996     else {
00997       RowMap = &(A().RowMatrixRowMap()); 
00998       ColMap = &(A().RowMatrixColMap()); 
00999       DomainMap = &(A().OperatorDomainMap());
01000     }
01001 
01002 #ifdef IFPACK_OVA_TIME_BUILD
01003     t1 = MPI_Wtime();
01004     fprintf(stderr,"[node %d]: overlap pointer pulls %2.3e\n", nodeID, t1-t0);
01005     t0=t1;
01006 #endif
01007     
01008     // For each column ID, determine the owning node (as opposed to core)
01009     // ID of the corresponding row.
01010     Epetra_IntVector colIdList( *ColMap );
01011     Epetra_IntVector rowIdList(*DomainMap);
01012     rowIdList.PutValue(nodeID);  
01013     Teuchos::RCP<Epetra_Import> nodeIdImporter = rcp(new Epetra_Import( *ColMap, *DomainMap ));
01014 
01015 #ifdef IFPACK_OVA_TIME_BUILD
01016     t1 = MPI_Wtime();
01017     fprintf(stderr,"[node %d]: overlap intvector/importer allocs %2.3e\n", nodeID, t1-t0);
01018     t0=t1;
01019 #endif
01020     
01021     colIdList.Import(rowIdList,*nodeIdImporter,Insert);
01022     
01023     int size = ColMap->NumMyElements() - RowMap->NumMyElements(); 
01024     int count = 0; 
01025 
01026 #ifdef IFPACK_OVA_TIME_BUILD
01027     t1 = MPI_Wtime();
01028     fprintf(stderr,"[node %d]: overlap col/row ID imports %2.3e\n", nodeID, t1-t0);
01029     t0=t1;
01030 #endif
01031 
01032     
01033     // define the set of off-node rows that are in ColMap but not in RowMap
01034     // This naturally does not include off-core rows that are on the same node as me, i.e., node buddy rows.
01035     for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) {
01036       int GID = ColMap->GID(i); 
01037       int myvotes=0;
01038       if (ghostTable.containsKey(GID)) myvotes = ghostTable.get(GID);
01039       if ( colIdList[i] != nodeID && myvotes < votesMax) // don't include anybody found in a previous round
01040       {
01041         int votes;
01042         if (ghostTable.containsKey(GID)) {
01043           votes = ghostTable.get(GID);
01044           votes++;
01045           ghostTable.put(GID,votes);
01046         } else {
01047           ghostTable.put(GID,1);
01048         }
01049       }
01050     } //for (int i = 0 ; i < ColMap->NumMyElements() ; ++i)
01051 
01052     Teuchos::Array<int> gidsarray,votesarray;
01053     ghostTable.arrayify(gidsarray,votesarray);
01054     int *gids = gidsarray.getRawPtr();
01055     int *votes = votesarray.getRawPtr();
01056 
01057     /*
01058        This next bit of code decides which node buddy (NB) gets which ghost points.  Everyone sends their
01059        list of ghost points to pid 0 of the local subcommunicator.  Pid 0 decides who gets what:
01060 
01061           - if a ghost point is touched by only one NB, that NB gets the ghost point
01062           - if two or more NBs share a ghost point, whichever NB has the most connections to the ghost
01063             point gets it.
01064     */
01065 
01066 #   ifdef HAVE_MPI  //FIXME What if we build in serial?!  This file will likely not compile.
01067     int *cullList;
01068     int ncull;
01069     //int mypid = nodeComm->MyPID();
01070 
01071     //FIXME timing
01072 #ifdef IFPACK_OVA_TIME_BUILD
01073     t1 = MPI_Wtime();
01074     fprintf(stderr,"[node %d]: overlap before culling %2.3e\n", nodeID, t1-t0);
01075     t0=t1;
01076 #endif
01077     //FIXME timing
01078 
01079 
01080     if (nodeComm->MyPID() == 0)
01081     {
01082       // Figure out how much pid 0 is to receive
01083       MPI_Status status;
01084       int *lengths = new int[nodeComm->NumProc()+1];
01085       lengths[0] = 0;
01086       lengths[1] = ghostTable.size();
01087       for (int i=1; i<nodeComm->NumProc(); i++) {
01088         int leng;
01089         MPI_Recv( &leng, 1, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
01090         lengths[i+1] = lengths[i] + leng;
01091       }
01092       int total = lengths[nodeComm->NumProc()];
01093 
01094       int* ghosts = new int[total];
01095       for (int i=0; i<total; i++) ghosts[i] = -9;
01096       int *round  = new int[total];
01097       int *owningpid  = new int[total];
01098 
01099       for (int i=1; i<nodeComm->NumProc(); i++) {
01100         int count = lengths[i+1] - lengths[i];
01101         MPI_Recv( ghosts+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
01102         MPI_Recv( round+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
01103       }
01104 
01105       // put in pid 0's info
01106       for (int i=0; i<lengths[1]; i++) {
01107         ghosts[i] = gids[i];
01108         round[i] = votes[i];
01109         owningpid[i] = 0;
01110       }
01111 
01112       // put in the pid associated with each ghost
01113       for (int j=1; j<nodeComm->NumProc(); j++) {
01114         for (int i=lengths[j]; i<lengths[j+1]; i++) {
01115           owningpid[i] = j;
01116         }
01117       }
01118 
01119       // sort everything based on the ghost gids
01120       int* roundpid[2];
01121       roundpid[0] = round; roundpid[1] = owningpid;
01122       Epetra_Util epetraUtil;
01123       epetraUtil.Sort(true,total,ghosts,0,0,2,roundpid);
01124 
01125       //set up arrays that get sent back to node buddies and that tell them what ghosts to cull
01126       int *nlosers = new int[nodeComm->NumProc()];
01127       int** losers = new int*[nodeComm->NumProc()];
01128       for (int i=0; i<nodeComm->NumProc(); i++) {
01129         nlosers[i] = 0;
01130         losers[i] = new int[ lengths[i+1]-lengths[i] ];
01131       }
01132 
01133       // Walk through ghosts array and and for each sequence of duplicate ghost GIDs, choose just one NB to keep it.
01134       // The logic is pretty simple.  The ghost list is sorted, so all duplicate PIDs are together.
01135       // The list is traversed.  As duplicates are found, node pid 0 keeps track of the current "winning"
01136       // pid.  When a pid is determined to have "lost" (less votes/connections to the current GID), the
01137       // GID is added to that pid's list of GIDs to be culled.  At the end of the repeated sequence, we have
01138       // a winner, and other NBs know whether they need to delete it from their import list.
01139       int max=0;   //for duplicated ghosts, index of pid with most votes and who hence keeps the ghost.
01140                    // TODO to break ties randomly
01141 
01142       for (int i=1; i<total; i++) {
01143         if (ghosts[i] == ghosts[i-1]) {
01144           int idx = i; // pid associated with idx is current "loser"
01145           if (round[i] > round[max]) {
01146             idx = max;
01147             max=i;
01148           }
01149           int j = owningpid[idx];
01150           losers[j][nlosers[j]++] = ghosts[idx];
01151         } else {
01152           max=i;
01153         }
01154       } //for (int i=1; i<total; i++)
01155 
01156       delete [] round;
01157       delete [] ghosts;
01158       delete [] owningpid;
01159 
01160       // send the arrays of ghost GIDs to be culled back to the respective node buddies
01161       for (int i=1; i<nodeComm->NumProc(); i++) {
01162         MPI_Send( nlosers+i, 1, MPI_INT, i, 8675, nodeComm->Comm());
01163         MPI_Send( losers[i], nlosers[i], MPI_INT, i, 8675, nodeComm->Comm());
01164       }
01165 
01166       //FIXME Unnecessary memory allocation and copying, but makes culling code cleaner
01167       //Could we stick this info into gids and votes, since neither is being used anymore?
01168       //TODO Instead of using "losers" arrays, just use "cullList" as in the else clause
01169       ncull = nlosers[0];
01170       cullList = new int[ncull+1];
01171       for (int i=0; i<nlosers[0]; i++)
01172         cullList[i] = losers[0][i];
01173 
01174       for (int i=0; i<nodeComm->NumProc(); i++)
01175         delete [] losers[i];
01176 
01177       delete [] lengths;
01178       delete [] losers;
01179       delete [] nlosers;
01180 
01181     } else { //everyone but pid 0
01182 
01183       // send to node pid 0 all ghosts that this pid could potentially import
01184       int hashsize = ghostTable.size();
01185       MPI_Send( &hashsize, 1, MPI_INT, 0, 8675, nodeComm->Comm());
01186       MPI_Send( gids, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
01187       MPI_Send( votes, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
01188 
01189       // receive the ghost GIDs that should not be imported (subset of the list sent off just above)
01190       MPI_Status status;
01191       MPI_Recv( &ncull, 1, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
01192       cullList = new int[ncull+1];
01193       MPI_Recv( cullList, ncull, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
01194     }
01195 
01196 
01197     //TODO clean out hash table after each time through overlap loop   4/1/07 JJH done moved both hash tables to local scope
01198 
01199     // Remove from my row map all off-node ghosts that will be imported by a node buddy.
01200     for (int i=0; i<ncull; i++) {
01201       try{ghostTable.remove(cullList[i]);}
01202 
01203       catch(...) {
01204         printf("pid %d: In OverlappingRowMatrix ctr, problem removing ghost elt %d from ghostTable\n",
01205                Comm().MyPID(),cullList[i]);
01206         fflush(stdout);
01207       }
01208     }//for
01209 
01210     delete [] cullList;
01211 
01212     // Save off the remaining ghost GIDs from the current overlap round.
01213     // These are off-node GIDs (rows) that I will import.
01214     gidsarray.clear(); votesarray.clear();
01215     ghostTable.arrayify(gidsarray,votesarray);
01216 
01217     std::vector<int> list(size); 
01218     count=0;
01219     for (int i=0; i<ghostTable.size(); i++) {
01220       // if votesarray[i] == votesmax, then this GID was found during a previous overlap round
01221       if (votesarray[i] < votesMax) {
01222         list[count] = gidsarray[i];
01223         ghostTable.put(gidsarray[i],votesMax); //set this GID's vote to the maximum so that this PID alway wins in any subsequent overlap tiebreaking.
01224         count++;
01225       }
01226     }
01227 
01228     //FIXME timing
01229 #ifdef IFPACK_OVA_TIME_BUILD
01230     t1 = MPI_Wtime();
01231     fprintf(stderr,"[node %d]: overlap duplicate removal %2.3e\n", nodeID, t1-t0);
01232     t0=t1;
01233 #endif
01234     //FIXME timing
01235 
01236 #   endif //ifdef HAVE_MPI
01237 
01238     TmpMap = rcp( new Epetra_Map(-1,count, &list[0],0,Comm()) );
01239 
01240     TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) ); 
01241 
01242     TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) ); 
01243 
01244     TmpMatrix->Import(A(),*TmpImporter,Insert); 
01245     TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap); 
01246 
01247     //FIXME timing
01248 #ifdef IFPACK_OVA_TIME_BUILD
01249     t1 = MPI_Wtime();
01250     fprintf(stderr,"[node %d]: overlap TmpMatrix fillcomplete %2.3e\n", nodeID, t1-t0);
01251     t0=t1;
01252 #endif
01253     //FIXME timing
01254 
01255 
01256     // These next two imports get the GIDs that need to go into the column map of the overlapped matrix.
01257 
01258     // For each column ID in the overlap, determine the owning node (as opposed to core)
01259     // ID of the corresponding row.  Save those column IDs whose owning node is the current one.
01260     // This should get all the imported ghost GIDs.
01261     Epetra_IntVector ov_colIdList( TmpMatrix->ColMap() );
01262     ov_colIdList.PutValue(-1);
01263     Epetra_IntVector ov_rowIdList( TmpMatrix->RowMap() );
01264     ov_rowIdList.PutValue(nodeID);  
01265     Teuchos::RCP<Epetra_Import> ov_nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), TmpMatrix->RowMap()));
01266     ov_colIdList.Import(ov_rowIdList,*ov_nodeIdImporter,Insert);
01267 
01268     for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
01269        if (ov_colIdList[i] == nodeID) {
01270          int GID = (ov_colIdList.Map()).GID(i);
01271          colMapTable.put(GID,1);
01272        }
01273     }
01274 
01275     // Do a second import of the owning node ID from A's rowmap to TmpMat's column map.  This ensures that
01276     // all GIDs that belong to a node buddy and are in a ghost row's sparsity pattern will be in the final
01277     // overlapped matrix's column map.
01278     ov_colIdList.PutValue(-1);
01279     Epetra_IntVector ArowIdList( A().RowMatrixRowMap() );
01280     ArowIdList.PutValue(nodeID);
01281     nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), A().RowMatrixRowMap() ));
01282     ov_colIdList.Import(ArowIdList,*nodeIdImporter,Insert);
01283 
01284     for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
01285        if (ov_colIdList[i] == nodeID) {
01286          int GID = (ov_colIdList.Map()).GID(i);
01287          colMapTable.put(GID,1);
01288        }
01289     }
01290 
01291     //FIXME timing
01292 #ifdef IFPACK_OVA_TIME_BUILD
01293     t1 = MPI_Wtime();
01294     fprintf(stderr,"[node %d]: overlap 2 imports to fix up colmap %2.3e\n", nodeID, t1-t0);
01295     t0=t1;
01296 #endif
01297     //FIXME timing
01298 
01299   } //for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap)
01300 
01301   /* ** ************************************************************************ ** */
01302   /* ** ********************** end of main overlap loop ************************ ** */
01303   /* ** ************************************************************************ ** */
01304 
01305   // off-node GIDs that will be in the overlap
01306   std::vector<int> ghostElements; 
01307 
01308   Teuchos::Array<int> gidsarray,votesarray;
01309   ghostTable.arrayify(gidsarray,votesarray);
01310   for (int i=0; i<ghostTable.size(); i++) {
01311     ghostElements.push_back(gidsarray[i]);
01312   }
01313 
01314     for (int i = 0 ; i < A().RowMatrixColMap().NumMyElements() ; ++i) {
01315       int GID = A().RowMatrixColMap().GID(i);
01316       // Remove any entries that are in A's original column map
01317       if (colMapTable.containsKey(GID)) {
01318         try{colMapTable.remove(GID);}
01319         catch(...) {
01320           printf("pid %d: In OverlappingRowMatrix ctr, problem removing entry %d from colMapTable\n", Comm().MyPID(),GID);
01321           fflush(stdout);
01322         }
01323       }
01324     }
01325 
01326     // GIDs that will be in the overlapped matrix's column map
01327     std::vector<int> colMapElements; 
01328 
01329   gidsarray.clear(); votesarray.clear();
01330   colMapTable.arrayify(gidsarray,votesarray);
01331   for (int i=0; i<colMapTable.size(); i++)
01332     colMapElements.push_back(gidsarray[i]);
01333 
01334 /*
01335    We need two maps here.  The first is the row map, which we've got by using my original rows
01336    plus whatever I've picked up in ghostElements.
01337 
01338    The second is a column map.  This map should include my rows, plus any columns that could come from node buddies.
01339    These GIDs come from the std:array colMapElements, which in turn comes from colMapTable.
01340    This map should *not* omit ghosts that have been imported by my node buddies, i.e., for any row that I own,
01341    the stencil should include all column GIDs (including imported ghosts) that are on the node.
01342 */
01343 
01344   // build the row map containing all the nodes (original matrix + off-node matrix)
01345   std::vector<int> rowList(NumMyRowsA_ + ghostElements.size());
01346   for (int i = 0 ; i < NumMyRowsA_ ; ++i)
01347     rowList[i] = A().RowMatrixRowMap().GID(i);
01348   for (int i = 0 ; i < (int)ghostElements.size() ; ++i)
01349     rowList[i + NumMyRowsA_] = ghostElements[i];
01350 
01351   // row map for the overlapped matrix (local + overlap)
01352   Map_ = rcp( new Epetra_Map(-1, NumMyRowsA_ + ghostElements.size(), &rowList[0], 0, Comm()) );
01353 
01354   // build the column map for the overlapping matrix
01355   //vector<int> colList(colMapElements.size());
01356   // column map for the overlapped matrix (local + overlap)
01357   //colMap_ = rcp( new Epetra_Map(-1, colMapElements.size(), &colList[0], 0, Comm()) );
01358   //for (int i = 0 ; i < (int)colMapElements.size() ; i++)
01359   //  colList[i] = colMapElements[i];
01360   std::vector<int> colList(A().RowMatrixColMap().NumMyElements() + colMapElements.size());
01361   int nc = A().RowMatrixColMap().NumMyElements();
01362   for (int i = 0 ; i < nc; i++)
01363     colList[i] = A().RowMatrixColMap().GID(i);
01364   for (int i = 0 ; i < (int)colMapElements.size() ; i++)
01365     colList[nc+i] = colMapElements[i];
01366 
01367   // column map for the overlapped matrix (local + overlap)
01368   //colMap_ = rcp( new Epetra_Map(-1, A().RowMatrixColMap().NumMyElements() + colMapElements.size(), &colList[0], 0, Comm()) );
01369   colMap_ = new Epetra_Map(-1, A().RowMatrixColMap().NumMyElements() + colMapElements.size(), &colList[0], 0, Comm());
01370     
01371 
01372   //FIXME timing
01373 #ifdef IFPACK_OVA_TIME_BUILD
01374   t1 = MPI_Wtime();
01375   fprintf(stderr,"[node %d]: time to init B() col/row maps %2.3e\n", nodeID, t1-t0);
01376   t0=t1;
01377 #endif
01378   //FIXME timing
01379 
01380   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
01381   //++++ start of sort
01382   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
01383   // build the column map, but don't use a copy constructor, b/c local communicator SubComm_ is
01384   // different from that of Matrix.
01385   try {
01386     // build row map based on the local communicator.  We need this temporarily to build the column map.
01387     Epetra_Map* nodeMap = new Epetra_Map(-1,NumMyRowsA_ + ghostElements.size(),&rowList[0],0,*nodeComm);
01388     int numMyElts = colMap_->NumMyElements();
01389     assert(numMyElts!=0);
01390 
01391     // The column map *must* be sorted: first locals, then ghosts.  The ghosts
01392     // must be further sorted so that they are contiguous by owning processor.
01393 
01394     // For each GID in column map, determine owning PID in local communicator.
01395     int* myGlobalElts = new int[numMyElts];
01396     colMap_->MyGlobalElements(myGlobalElts);
01397     int* pidList = new int[numMyElts];
01398     nodeMap->RemoteIDList(numMyElts, myGlobalElts, pidList, 0);
01399 
01400     // First sort column map GIDs based on the owning PID in local communicator.
01401     Epetra_Util Util;
01402     int *tt[1];
01403     tt[0] = myGlobalElts;
01404     Util.Sort(true, numMyElts, pidList, 0, (double**)0, 1, tt);
01405 
01406     // For each remote PID, sort the column map GIDs in ascending order.
01407     // Don't sort the local PID's entries.
01408     int localStart=0;
01409     while (localStart<numMyElts) {
01410       int currPID = (pidList)[localStart];
01411       int i=localStart;
01412       while (i<numMyElts && (pidList)[i] == currPID) i++;
01413       if (currPID != nodeComm->MyPID())
01414         Util.Sort(true, i-localStart, (myGlobalElts)+localStart, 0, 0, 0, 0);
01415       localStart = i;
01416     }
01417 
01418     // now move the local entries to the front of the list
01419     localStart=0;
01420     while (localStart<numMyElts && (pidList)[localStart] != nodeComm->MyPID()) localStart++;
01421     assert(localStart != numMyElts);
01422     int localEnd=localStart;
01423     while (localEnd<numMyElts && (pidList)[localEnd] == nodeComm->MyPID()) localEnd++;
01424     int* mySortedGlobalElts = new int[numMyElts];
01425     int leng = localEnd - localStart;
01426     /* This is a little gotcha.  It appears that the ordering of the column map's local entries
01427        must be the same as that of the domain map's local entries.  See the comment in method
01428        MakeColMap() in Epetra_CrsGraph.cpp, line 1072. */
01429     int *rowGlobalElts =  nodeMap->MyGlobalElements();
01430 
01431     /*TODO TODO TODO NTS rows 68 and 83 show up as local GIDs in rowGlobalElts for both pids 0 & 1 on node 0.
01432                     This seems to imply that there is something wrong with rowList!!
01433                     It is ok for the overlapped matrix's row map to have repeated entries (it's overlapped, after all).
01434                     But ... on a node, there should be no repeats.
01435                     Ok, here's the underlying cause.  On node 0, gpid 1 finds 83 on overlap round 0.  gpid 0 finds 83
01436                     on overlap round 1.  This shouldn't happen .. once a node buddy "owns" a row, no other node buddy
01437                     should ever import ("own") that row.
01438 
01439                     Here's a possible fix.  Extend tie breaking across overlap rounds.  This means sending to lpid 0
01440                     a list of *all* rows you've imported (this round plus previous rounds) for tie-breaking
01441                     If a nb won a ghost row in a previous round, his votes for that ghost row should be really high, i.e.,
01442                     he should probably win that ghost row forever.
01443                     7/17/09
01444 
01445       7/28/09 JJH I believe I've fixed this, but I thought it might be helpful to have this comment in here, in case I missed something.
01446     */
01447 
01448     //move locals to front of list
01449     for (int i=0; i<leng; i++) {
01450       /* //original code */
01451       (mySortedGlobalElts)[i] = rowGlobalElts[i];
01452       //(&*mySortedGlobalElts)[i] = (&*myGlobalElts)[localStart+i];
01453       //(&*mySortedPidList)[i] = (&*pidList)[localStart+i];
01454     }
01455     for (int i=leng; i< localEnd; i++) {
01456       (mySortedGlobalElts)[i] = (myGlobalElts)[i-leng];
01457     }
01458     for (int i=localEnd; i<numMyElts; i++) {
01459       (mySortedGlobalElts)[i] = (myGlobalElts)[i];
01460     }
01461 
01462     //FIXME timing
01463 #ifdef IFPACK_OVA_TIME_BUILD
01464     t1 = MPI_Wtime();
01465     fprintf(stderr,"[node %d]: time to sort col map arrays (cp=1) %2.3e\n", nodeID, t1-t0);
01466     t0=t1;
01467 #endif
01468     //FIXME timing
01469 
01470     int indexBase = colMap_->IndexBase();
01471     if (colMap_) delete colMap_;
01472     colMap_ = new Epetra_Map(-1,numMyElts,mySortedGlobalElts,indexBase,Comm());
01473 
01474     delete nodeMap;
01475     delete [] myGlobalElts;
01476     delete [] pidList;
01477     delete [] mySortedGlobalElts;
01478 
01479   } //try
01480   catch(...) {
01481     printf("** * Ifpack_OverlappingRowmatrix ctor: problem creating column map * **\n\n");
01482   }
01483 
01484   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
01485   //++++ end of sort
01486   //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
01487 
01488   //FIXME timing
01489 #ifdef IFPACK_OVA_TIME_BUILD
01490   t1 = MPI_Wtime();
01491   fprintf(stderr,"[node %d]: time to sort col map arrays (cp=2) %2.3e\n", nodeID, t1-t0);
01492   t0=t1;
01493 #endif
01494   //FIXME timing
01495 
01496 
01497 /*
01498 
01499    FIXME
01500    Does the column map need to be sorted for the overlapping matrix?
01501 
01502    The column map *must* be sorted:
01503 
01504         first locals
01505         then ghosts
01506 
01507    The ghosts must be further sorted so that they are contiguous by owning processor
01508 
01509   int* RemoteSizeList = 0
01510   int* RemoteColIndices = ColIndices.Values() + NumLocalColGIDs; // Points to back end of ColIndices
01511 
01512   EPETRA_CHK_ERR(DomainMap.RemoteIDList(NumRemoteColGIDs, RemoteColIndices, PIDList.Values(), 0, RemoteSizeList));
01513   Epetra_Util epetraUtil;
01514   SortLists[0] = RemoteColIndices;
01515   SortLists[1] = RemoteSizeList;
01516   epetraUtil.Sort(true, NumRemoteColGIDs, PIDList.Values(), 0, 0, NLists, SortLists);
01517 */
01518 
01519   // now build the map corresponding to all the external nodes
01520   // (with respect to A().RowMatrixRowMap().
01521   ExtMap_ = rcp( new Epetra_Map(-1,ghostElements.size(), &ghostElements[0],0,Comm()) );
01522   ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*colMap_,0) ); 
01523 
01524   ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) ); 
01525   ExtMatrix_->Import(A(),*ExtImporter_,Insert); 
01526 
01527   //FIXME timing
01528 #ifdef IFPACK_OVA_TIME_BUILD
01529   t1 = MPI_Wtime();
01530   fprintf(stderr,"[node %d]: time to import and setup ExtMap_ %2.3e\n", nodeID, t1-t0);
01531   t0=t1;
01532 #endif
01533   //FIXME timing
01534 
01535 
01536 /*
01537   Notes to self:    In FillComplete, the range map does not have to be 1-1 as long as
01538                     (row map == range map).  Ditto for the domain map not being 1-1
01539                     if (col map == domain map).
01540                     
01541 */
01542 
01543   ExtMatrix_->FillComplete( *colMap_ , *Map_ ); //FIXME wrong
01544 
01545   //FIXME timing
01546 #ifdef IFPACK_OVA_TIME_BUILD
01547   t1 = MPI_Wtime();
01548   fprintf(stderr,"[node %d]: time to FillComplete B() %2.3e\n", nodeID, t1-t0);
01549   t0=t1;
01550 #endif
01551   //FIXME timing
01552 
01553 
01554   // Note: B() = *ExtMatrix_ .
01555 
01556   Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) ); //FIXME is this right?!
01557 
01558   // fix indices for overlapping matrix
01559   NumMyRowsB_ = B().NumMyRows();
01560   NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;  //TODO is this wrong for a subdomain on >1 processor? // should be ok
01561   //NumMyCols_ = NumMyRows_;  //TODO is this wrong for a subdomain on >1 processor?  // YES!!!
01562   //NumMyCols_ = A().NumMyCols() + B().NumMyCols();
01563   NumMyCols_ = B().NumMyCols();
01564 
01565   /*FIXME*/ //somehow B's NumMyCols is the entire subdomain (local + overlap)
01566 
01567   NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
01568   
01569   NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
01570   Comm().SumAll(&NumMyNonzeros_,&NumGlobalNonzeros_,1);
01571   MaxNumEntries_ = A().MaxNumEntries();
01572   
01573   if (MaxNumEntries_ < B().MaxNumEntries())
01574     MaxNumEntries_ = B().MaxNumEntries();
01575 
01576 # ifdef HAVE_MPI
01577   delete nodeComm;
01578 # endif
01579 
01580   //FIXME timing
01581 #ifdef IFPACK_OVA_TIME_BUILD
01582   t1 = MPI_Wtime();
01583   fprintf(stderr,"[node %d]: time for final calcs %2.3e\n", nodeID, t1-t0);
01584   fprintf(stderr,"[node %d]: total IORM ctor time %2.3e\n", nodeID, t1-true_t0);
01585 #endif
01586   //FIXME timing
01587 
01588 
01589 } //Ifpack_OverlappingRowMatrix() ctor for more than one core
01590 
01591 // Destructor
01592 Ifpack_OverlappingRowMatrix::~Ifpack_OverlappingRowMatrix() {
01593   delete colMap_;
01594 }
01595 
01596 # endif //ifdef IFPACK_NODE_AWARE_CODE
01597 #endif // ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01598 
01599 
01600 // ======================================================================
01601 int Ifpack_OverlappingRowMatrix::
01602 NumMyRowEntries(int MyRow, int & NumEntries) const
01603 {
01604   if (MyRow < NumMyRowsA_)
01605     return(A().NumMyRowEntries(MyRow,NumEntries));
01606   else
01607     return(B().NumMyRowEntries(MyRow - NumMyRowsA_, NumEntries));
01608 }
01609 
01610 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01611 // ======================================================================
01612 int Ifpack_OverlappingRowMatrix::
01613 ExtractMyRowCopy(int LocRow, int Length, int & NumEntries, double *Values, 
01614                  int * Indices) const
01615 {
01616   //assert(1==0);
01617   int ierr;
01618   const Epetra_Map *Themap;
01619   if (LocRow < NumMyRowsA_) {
01620     ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01621     Themap=&A().RowMatrixColMap();
01622   }
01623   else {
01624     ierr = B().ExtractMyRowCopy(LocRow-NumMyRowsA_,Length,NumEntries,Values,Indices);
01625     Themap=&B().RowMatrixColMap();
01626   }
01627 
01628   IFPACK_RETURN(ierr);
01629 }
01630 
01631 // ======================================================================
01632 int Ifpack_OverlappingRowMatrix::
01633 ExtractGlobalRowCopy(int GlobRow, int Length, int & NumEntries, double *Values, 
01634                  int * Indices) const
01635 {
01636   int ierr;
01637   const Epetra_Map *Themap;
01638   int LocRow = A().RowMatrixRowMap().LID(GlobRow);
01639   if (LocRow < NumMyRowsA_ && LocRow != -1) { //TODO don't need to check less than nummyrows
01640     ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01641     Themap=&A().RowMatrixColMap();
01642   }
01643   else {
01644     LocRow = B().RowMatrixRowMap().LID(GlobRow);
01645     assert(LocRow!=-1);
01646     //ierr = B().ExtractMyRowCopy(LocRow-NumMyRowsA_,Length,NumEntries,Values,Indices);
01647     ierr = B().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01648     Themap=&B().RowMatrixColMap();
01649   }
01650 
01651   for (int i=0; i<NumEntries; i++) {
01652     Indices[i]=Themap->GID(Indices[i]);
01653     assert(Indices[i] != -1);
01654   }
01655 
01656   IFPACK_RETURN(ierr);
01657 }
01658 #else
01659 # ifdef IFPACK_NODE_AWARE_CODE
01660 // ======================================================================
01661 int Ifpack_OverlappingRowMatrix::
01662 ExtractMyRowCopy(int LocRow, int Length, int & NumEntries, double *Values, 
01663                  int * Indices) const
01664 {
01665   assert(1==0);
01666   int ierr;
01667   const Epetra_Map *Themap;
01668   if (LocRow < NumMyRowsA_) {
01669     ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01670     Themap=&A().RowMatrixColMap();
01671   }
01672   else {
01673     ierr = B().ExtractMyRowCopy(LocRow-NumMyRowsA_,Length,NumEntries,Values,Indices);
01674     Themap=&B().RowMatrixColMap();
01675   }
01676 
01677   IFPACK_RETURN(ierr);
01678 }
01679 
01680 // ======================================================================
01681 int Ifpack_OverlappingRowMatrix::
01682 ExtractGlobalRowCopy(int GlobRow, int Length, int & NumEntries, double *Values, 
01683                  int * Indices) const
01684 {
01685   int ierr;
01686   const Epetra_Map *Themap;
01687   int LocRow = A().RowMatrixRowMap().LID(GlobRow);
01688   if (LocRow < NumMyRowsA_ && LocRow != -1) { //TODO don't need to check less than nummyrows
01689     ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01690     Themap=&A().RowMatrixColMap();
01691   }
01692   else {
01693     LocRow = B().RowMatrixRowMap().LID(GlobRow);
01694     assert(LocRow!=-1);
01695     //ierr = B().ExtractMyRowCopy(LocRow-NumMyRowsA_,Length,NumEntries,Values,Indices);
01696     ierr = B().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01697     Themap=&B().RowMatrixColMap();
01698   }
01699 
01700   for (int i=0; i<NumEntries; i++) {
01701     Indices[i]=Themap->GID(Indices[i]);
01702     assert(Indices[i] != -1);
01703   }
01704 
01705   IFPACK_RETURN(ierr);
01706 }
01707 # else
01708 
01709 // ======================================================================
01710 int Ifpack_OverlappingRowMatrix::
01711 ExtractMyRowCopy(int MyRow, int Length, int & NumEntries, double *Values,
01712                  int * Indices) const
01713 {
01714   int ierr;
01715   if (MyRow < NumMyRowsA_)
01716     ierr = A().ExtractMyRowCopy(MyRow,Length,NumEntries,Values,Indices);
01717   else
01718     ierr = B().ExtractMyRowCopy(MyRow - NumMyRowsA_,Length,NumEntries,
01719                                 Values,Indices);
01720 
01721   IFPACK_RETURN(ierr);
01722 }
01723 # endif // ifdef IFPACK_NODE_AWARE_CODE
01724 #endif // ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01725 
01726 // ======================================================================
01727 int Ifpack_OverlappingRowMatrix::
01728 ExtractDiagonalCopy(Epetra_Vector & Diagonal) const
01729 {
01730   IFPACK_CHK_ERR(-1);
01731 }
01732 
01733 
01734 // ======================================================================
01735 int Ifpack_OverlappingRowMatrix::
01736 Multiply(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
01737 {
01738   int NumVectors = X.NumVectors();
01739   std::vector<int> Ind(MaxNumEntries_);
01740   std::vector<double> Val(MaxNumEntries_);
01741 
01742   Y.PutScalar(0.0);
01743 
01744   // matvec with A (local rows)
01745   for (int i = 0 ; i < NumMyRowsA_ ; ++i) {
01746     for (int k = 0 ; k < NumVectors ; ++k) {
01747       int Nnz;
01748       IFPACK_CHK_ERR(A().ExtractMyRowCopy(i,MaxNumEntries_,Nnz, 
01749                                           &Val[0], &Ind[0]));
01750       for (int j = 0 ; j < Nnz ; ++j) {
01751         Y[k][i] += Val[j] * X[k][Ind[j]];
01752       }
01753     }
01754   }
01755 
01756   // matvec with B (overlapping rows)
01757   for (int i = 0 ; i < NumMyRowsB_ ; ++i) {
01758     for (int k = 0 ; k < NumVectors ; ++k) {
01759       int Nnz;
01760       IFPACK_CHK_ERR(B().ExtractMyRowCopy(i,MaxNumEntries_,Nnz, 
01761                                           &Val[0], &Ind[0]));
01762       for (int j = 0 ; j < Nnz ; ++j) {
01763         Y[k][i + NumMyRowsA_] += Val[j] * X[k][Ind[j]];
01764       }
01765     }
01766   }
01767   return(0);
01768 }
01769 
01770 // ======================================================================
01771 int Ifpack_OverlappingRowMatrix::
01772 Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
01773 {
01774   IFPACK_CHK_ERR(Multiply(UseTranspose(),X,Y));
01775   return(0);
01776 }
01777 
01778 // ======================================================================
01779 int Ifpack_OverlappingRowMatrix::
01780 ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
01781 {
01782   IFPACK_CHK_ERR(-1);
01783 }
01784 
01785 // ======================================================================
01786 #ifndef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01787 # ifndef IFPACK_NODE_AWARE_CODE
01788 Epetra_RowMatrix& Ifpack_OverlappingRowMatrix::B() const
01789 {
01790   return(*ExtMatrix_);
01791 }
01792 # endif
01793 #endif
01794 // ======================================================================
01795 const Epetra_BlockMap& Ifpack_OverlappingRowMatrix::Map() const
01796 {
01797   return(*Map_);
01798 }
01799 
01800 // ======================================================================
01801 int Ifpack_OverlappingRowMatrix::
01802 ImportMultiVector(const Epetra_MultiVector& X, Epetra_MultiVector& OvX,
01803                   Epetra_CombineMode CM)
01804 {
01805   OvX.Import(X,*Importer_,CM);
01806   return(0);
01807 }
01808 
01809 // ======================================================================
01810 int Ifpack_OverlappingRowMatrix::
01811 ExportMultiVector(const Epetra_MultiVector& OvX, Epetra_MultiVector& X,
01812                   Epetra_CombineMode CM)
01813 {
01814   X.Export(OvX,*Importer_,CM);
01815   return(0);
01816 }
01817 
 All Classes Files Functions Variables Enumerations Friends