//============================================================================== // Epetra_Export constructor for a Epetra_BlockMap object Epetra_Export::Epetra_Export( const Epetra_BlockMap & SourceMap, const Epetra_BlockMap & TargetMap) : Epetra_Object("Epetra::Export"), TargetMap_(TargetMap), SourceMap_(SourceMap), NumSameIDs_(0), NumPermuteIDs_(0), PermuteToLIDs_(0), PermuteFromLIDs_(0), NumRemoteIDs_(0), RemoteLIDs_(0), NumExportIDs_(0), ExportLIDs_(0), ExportPIDs_(0), NumSend_(0), NumRecv_(0), Distor_(0) { int i; // Build three ID lists: // NumSameIDs - Number of IDs in TargetMap and SourceMap that are identical, up to the first // nonidentical ID. // NumPermuteIDs - Number of IDs in SourceMap that must be indirectly loaded but are on this processor. // NumExportIDs - Number of IDs that are in SourceMap but not in TargetMap, and thus must be exported. int NumSourceIDs = SourceMap.NumMyElements(); int NumTargetIDs = TargetMap.NumMyElements(); int *TargetGIDs = 0; if (NumTargetIDs>0) { TargetGIDs = new int[NumTargetIDs]; TargetMap.MyGlobalElements(TargetGIDs); } int * SourceGIDs = 0; if (NumSourceIDs>0) { SourceGIDs = new int[NumSourceIDs]; SourceMap.MyGlobalElements(SourceGIDs); } int MinIDs = EPETRA_MIN(NumSourceIDs, NumTargetIDs); NumSameIDs_ = 0; for (i=0; i< MinIDs; i++) if (TargetGIDs[i]==SourceGIDs[i]) NumSameIDs_++; else break; // Find count of Source IDs that are truly remote and those that are local but permuted NumPermuteIDs_ = 0; NumExportIDs_ = 0; for (i=NumSameIDs_; i< NumSourceIDs; i++) if (TargetMap.MyGID(SourceGIDs[i])) NumPermuteIDs_++; // Check if Source GID is a local Target GID else NumExportIDs_++; // If not, then it is remote // Define remote and permutation lists int * ExportGIDs = 0; if (NumExportIDs_>0) { ExportLIDs_ = new int[NumExportIDs_]; ExportGIDs = new int[NumExportIDs_]; } if (NumPermuteIDs_>0) { PermuteToLIDs_ = new int[NumPermuteIDs_]; PermuteFromLIDs_ = new int[NumPermuteIDs_]; } NumPermuteIDs_ = 0; NumExportIDs_ = 0; for (i=NumSameIDs_; i< NumSourceIDs; i++) { if (TargetMap.MyGID(SourceGIDs[i])) { PermuteFromLIDs_[NumPermuteIDs_] = i; PermuteToLIDs_[NumPermuteIDs_++] = TargetMap.LID(SourceGIDs[i]); } else { //NumSend_ +=SourceMap.ElementSize(i); // Count total number of entries to send NumSend_ +=SourceMap.MaxElementSize(); // Count total number of entries to send (currently need max) ExportGIDs[NumExportIDs_] = SourceGIDs[i]; ExportLIDs_[NumExportIDs_++] = i; } } if ( NumExportIDs_>0 && !SourceMap.DistributedGlobal()) ReportError("Warning in Epetra_Export: Serial Export has remote IDs. (Exporting from Subset of Source Map)", 1); // Test for distributed cases int ierr = 0; if (SourceMap.DistributedGlobal()) { if (NumExportIDs_>0) ExportPIDs_ = new int[NumExportIDs_]; ierr = TargetMap.RemoteIDList(NumExportIDs_, ExportGIDs, ExportPIDs_, 0); // Get remote PIDs if( ierr ) throw ReportError("Error in Epetra_BlockMap::RemoteIDList", ierr); //Get rid of IDs not in Target Map if(NumExportIDs_>0) { int cnt = 0; for( i = 0; i < NumExportIDs_; ++i ) if( ExportPIDs_[i] == -1 ) ++cnt; if( cnt ) { int * NewExportGIDs = 0; int * NewExportPIDs = 0; int * NewExportLIDs = 0; int cnt1 = NumExportIDs_-cnt; if (cnt1) { NewExportGIDs = new int[cnt1]; NewExportPIDs = new int[cnt1]; NewExportLIDs = new int[cnt1]; } cnt = 0; for( i = 0; i < NumExportIDs_; ++i ) if( ExportPIDs_[i] != -1 ) { NewExportGIDs[cnt] = ExportGIDs[i]; NewExportPIDs[cnt] = ExportPIDs_[i]; NewExportLIDs[cnt] = ExportLIDs_[i]; ++cnt; } assert(cnt==cnt1); // Sanity test NumExportIDs_ = cnt; delete [] ExportGIDs; delete [] ExportPIDs_; delete [] ExportLIDs_; ExportGIDs = NewExportGIDs; ExportPIDs_ = NewExportPIDs; ExportLIDs_ = NewExportLIDs; ReportError("Warning in Epetra_Export: Source IDs not found in Target Map (Do you want to export from subset of Source Map?)", 1 ); } } //Make sure Export IDs are ordered by processor Epetra_Util util; int * tmpPtr[2]; tmpPtr[0] = ExportLIDs_, tmpPtr[1] = ExportGIDs; util.Sort(true,NumExportIDs_,ExportPIDs_,0,0,2,tmpPtr); Distor_ = SourceMap.Comm().CreateDistributor(); // Construct list of exports that calling processor needs to send as a result // of everyone asking for what it needs to receive. ierr = Distor_->CreateFromSends( NumExportIDs_, ExportPIDs_, true, NumRemoteIDs_); if (ierr!=0) throw ReportError("Error in Epetra_Distributor.CreateFromSends()", ierr); // Use comm plan with ExportGIDs to find out who is sending to us and // get proper ordering of GIDs for remote entries // (that we will convert to LIDs when done). if (NumRemoteIDs_>0) RemoteLIDs_ = new int[NumRemoteIDs_]; // Allocate space for LIDs in target that are // going to get something from off-processor. char * cRemoteGIDs = 0; //Do will alloc memory for this object int LenCRemoteGIDs = 0; ierr = Distor_->Do(reinterpret_cast<char *> (ExportGIDs), sizeof( int ), LenCRemoteGIDs, cRemoteGIDs); if (ierr) throw ReportError("Error in Epetra_Distributor.Do()", ierr); int * RemoteGIDs = reinterpret_cast<int*>(cRemoteGIDs); // Remote IDs come in as GIDs, convert to LIDs for (i=0; i< NumRemoteIDs_; i++) { RemoteLIDs_[i] = TargetMap.LID(RemoteGIDs[i]); //NumRecv_ += TargetMap.ElementSize(RemoteLIDs_[i]); // Count total number of entries to receive NumRecv_ += TargetMap.MaxElementSize(); // Count total number of entries to receive (currently need max) } if (NumExportIDs_>0) delete [] ExportGIDs; if (LenCRemoteGIDs>0) delete [] cRemoteGIDs; } if (NumTargetIDs>0) delete [] TargetGIDs; if (NumSourceIDs>0) delete [] SourceGIDs; return; }
int checkmap(Epetra_BlockMap & Map, int NumGlobalElements, int NumMyElements, int *MyGlobalElements, int ElementSize, int * ElementSizeList, int NumGlobalPoints, int NumMyPoints, int IndexBase, Epetra_Comm& Comm, bool DistributedGlobal, bool IsOneToOne) { int i, ierr=0, forierr=0;// forierr is used in for loops, then is tested // after for loop completes to see if it is non zero - potentially prevents // thousands of error messages if (ElementSizeList==0) { EPETRA_TEST_ERR(!Map.ConstantElementSize(),ierr); } else EPETRA_TEST_ERR(Map.ConstantElementSize(),ierr); EPETRA_TEST_ERR(DistributedGlobal!=Map.DistributedGlobal(),ierr); EPETRA_TEST_ERR(IsOneToOne!=Map.IsOneToOne(),ierr); int *MyElementSizeList; if (ElementSizeList==0) { EPETRA_TEST_ERR(Map.ElementSize()!=ElementSize,ierr); MyElementSizeList = new int[NumMyElements]; EPETRA_TEST_ERR(Map.ElementSizeList(MyElementSizeList)!=0,ierr); forierr = 0; for (i=0; i<NumMyElements; i++) forierr += MyElementSizeList[i]!=ElementSize; EPETRA_TEST_ERR(forierr,ierr); EPETRA_TEST_ERR(Map.MaxMyElementSize() != ElementSize,ierr); EPETRA_TEST_ERR(Map.MinMyElementSize() != ElementSize,ierr); } else { MyElementSizeList = new int[NumMyElements]; EPETRA_TEST_ERR(Map.ElementSizeList(MyElementSizeList)!=0,ierr); int MaxSize = MyElementSizeList[0]; int MinSize = MyElementSizeList[0]; forierr=0; for (i=0; i<NumMyElements; i++) { forierr += MyElementSizeList[i]!=ElementSizeList[i]; if (MyElementSizeList[i] > MaxSize) MaxSize = MyElementSizeList[i]; if (MyElementSizeList[i] < MinSize) MinSize = MyElementSizeList[i]; // Test ElementSize(int LID) method forierr += Map.ElementSize(Map.LID(MyGlobalElements[i])) != ElementSizeList[i]; } EPETRA_TEST_ERR(forierr,ierr); EPETRA_TEST_ERR(MaxSize !=Map.MaxMyElementSize(),ierr); EPETRA_TEST_ERR(MinSize !=Map.MinMyElementSize(),ierr); } const Epetra_Comm & Comm1 = Map.Comm(); EPETRA_TEST_ERR(Comm1.NumProc()!=Comm.NumProc(),ierr); EPETRA_TEST_ERR(Comm1.MyPID()!=Comm.MyPID(),ierr); EPETRA_TEST_ERR(Map.IndexBase()!=IndexBase,ierr); EPETRA_TEST_ERR(!Map.LinearMap() && MyGlobalElements==0,ierr); EPETRA_TEST_ERR(Map.LinearMap() && MyGlobalElements!=0,ierr); EPETRA_TEST_ERR(Map.MaxAllGID()!=NumGlobalElements-1+IndexBase,ierr); EPETRA_TEST_ERR(Map.MaxElementSize()!=ElementSize,ierr); int MaxLID = Map.MaxLID(); EPETRA_TEST_ERR(MaxLID!=NumMyElements-1,ierr); int MaxMyGID = (Comm.MyPID()+1)*NumMyElements-1+IndexBase; if (Comm.MyPID()>2) MaxMyGID+=3; if (!DistributedGlobal) MaxMyGID = NumMyElements-1+IndexBase; EPETRA_TEST_ERR(Map.MaxMyGID()!=MaxMyGID,ierr); EPETRA_TEST_ERR(Map.MinAllGID()!=IndexBase,ierr); if (ElementSizeList==0) { EPETRA_TEST_ERR(Map.MinElementSize()!=ElementSize,ierr); } else EPETRA_TEST_ERR(Map.MinElementSize()!=2,ierr); int MinLID = Map.MinLID(); EPETRA_TEST_ERR(MinLID!=0,ierr); int MinMyGID = Comm.MyPID()*NumMyElements+IndexBase; if (Comm.MyPID()>2) MinMyGID+=3; if (!DistributedGlobal) MinMyGID = IndexBase; // Not really needed EPETRA_TEST_ERR(Map.MinMyGID()!=MinMyGID,ierr); int * MyGlobalElements1 = new int[NumMyElements]; EPETRA_TEST_ERR(Map.MyGlobalElements(MyGlobalElements1)!=0,ierr); forierr = 0; if (MyGlobalElements==0) { for (i=0; i<NumMyElements; i++) forierr += MyGlobalElements1[i]!=MinMyGID+i; EPETRA_TEST_ERR(forierr,ierr); } else { for (i=0; i<NumMyElements; i++) forierr += MyGlobalElements[i]!=MyGlobalElements1[i]; EPETRA_TEST_ERR(forierr,ierr); } EPETRA_TEST_ERR(Map.NumGlobalElements()!=NumGlobalElements,ierr); EPETRA_TEST_ERR(Map.NumGlobalPoints()!=NumGlobalPoints,ierr); EPETRA_TEST_ERR(Map.NumMyElements()!=NumMyElements,ierr); EPETRA_TEST_ERR(Map.NumMyPoints()!=NumMyPoints,ierr); int MaxMyGID2 = Map.GID(Map.LID(MaxMyGID)); EPETRA_TEST_ERR(MaxMyGID2 != MaxMyGID,ierr); int MaxLID2 = Map.LID(Map.GID(MaxLID)); EPETRA_TEST_ERR(MaxLID2 != MaxLID,ierr); EPETRA_TEST_ERR(Map.GID(MaxLID+1) != IndexBase-1,ierr);// MaxLID+1 doesn't exist EPETRA_TEST_ERR(Map.LID(MaxMyGID+1) != -1,ierr);// MaxMyGID+1 doesn't exist or is on a different processor EPETRA_TEST_ERR(!Map.MyGID(MaxMyGID),ierr); EPETRA_TEST_ERR(Map.MyGID(MaxMyGID+1),ierr); EPETRA_TEST_ERR(!Map.MyLID(MaxLID),ierr); EPETRA_TEST_ERR(Map.MyLID(MaxLID+1),ierr); EPETRA_TEST_ERR(!Map.MyGID(Map.GID(MaxLID)),ierr); EPETRA_TEST_ERR(Map.MyGID(Map.GID(MaxLID+1)),ierr); EPETRA_TEST_ERR(!Map.MyLID(Map.LID(MaxMyGID)),ierr); EPETRA_TEST_ERR(Map.MyLID(Map.LID(MaxMyGID+1)),ierr); // Test the FirstPointInElementList methods, begin by testing that they produce identical results int * FirstPointInElementList = new int[NumMyElements+1]; Map.FirstPointInElementList(FirstPointInElementList); int * FirstPointInElementList1 = Map.FirstPointInElementList(); forierr = 0; for (i=0; i<=NumMyElements; i++) forierr += FirstPointInElementList[i]!=FirstPointInElementList1[i]; EPETRA_TEST_ERR(forierr,ierr); // Now make sure values are correct forierr = 0; if (Map.ConstantElementSize()) { for (i=0; i<=NumMyElements; i++) forierr += FirstPointInElementList1[i]!=(i*ElementSize);// NOTE:FirstPointInElement[NumMyElements] is not the first point of an element EPETRA_TEST_ERR(forierr,ierr); } else { int FirstPoint = 0; for (i=0; i<NumMyElements; i++) { forierr += FirstPointInElementList1[i]!=FirstPoint; FirstPoint += ElementSizeList[i]; } EPETRA_TEST_ERR(forierr,ierr); EPETRA_TEST_ERR(FirstPointInElementList[NumMyElements] != NumMyPoints,ierr);// The last entry in the array = the total number of Points on the proc } delete [] FirstPointInElementList; // Declare some variables for the FindLocalElementID test int ElementID, Offset; // Test the PointToElementList methods, begin by testing that they produce identical results int * PointToElementList = new int[NumMyPoints]; Map.PointToElementList(PointToElementList); int * PointToElementList1 = Map.PointToElementList(); forierr = 0; for (i=0; i<NumMyPoints; i++) forierr += PointToElementList1[i] != PointToElementList[i]; EPETRA_TEST_ERR(forierr,ierr); //Now make sure values are correct forierr=0; if (Map.ConstantElementSize()) { for (i=0; i<NumMyElements; i++) for (int j=0; j<ElementSize; j++) { forierr += PointToElementList[i*ElementSize+j] != i; // Test FindLocalElementID method Map.FindLocalElementID(i*ElementSize+j,ElementID,Offset); forierr += ElementID != i || Offset != j; } EPETRA_TEST_ERR(forierr,ierr); } else { int MyPointTot = 0; // Keep track of total number of points in all previously completely checked elements for (i=0; i<NumMyElements; i++) { for (int j=0; j<ElementSizeList[i]; j++) { forierr += PointToElementList[MyPointTot+j] != i; // Test FindLocalElementID method Map.FindLocalElementID(MyPointTot+j,ElementID,Offset); forierr += ElementID != i || Offset != j; } MyPointTot += ElementSizeList[i]; } EPETRA_TEST_ERR(forierr,ierr); } delete [] PointToElementList; // Check RemoteIDList function that includes a parameter for size // Get some GIDs off of each processor to test int TotalNumEle, NumElePerProc, NumProc = Comm.NumProc(); int MinNumEleOnProc; int NumMyEle = Map.NumMyElements(); Comm.MinAll(&NumMyEle,&MinNumEleOnProc,1); if (MinNumEleOnProc > 5) NumElePerProc = 6; else NumElePerProc = MinNumEleOnProc; if (NumElePerProc > 0) { TotalNumEle = NumElePerProc*NumProc; int * MyGIDlist = new int[NumElePerProc]; int * GIDlist = new int[TotalNumEle]; int * PIDlist = new int[TotalNumEle]; int * LIDlist = new int[TotalNumEle]; int * SizeList = new int[TotalNumEle]; for (i=0; i<NumElePerProc; i++) MyGIDlist[i] = MyGlobalElements1[i]; Comm.GatherAll(MyGIDlist,GIDlist,NumElePerProc);// Get a few values from each proc Map.RemoteIDList(TotalNumEle, GIDlist, PIDlist, LIDlist, SizeList); int MyPID= Comm.MyPID(); forierr = 0; for (i=0; i<TotalNumEle; i++) { if (Map.MyGID(GIDlist[i])) { forierr += PIDlist[i] != MyPID; forierr += !Map.MyLID(Map.LID(GIDlist[i])) || Map.LID(GIDlist[i]) != LIDlist[i] || Map.GID(LIDlist[i]) != GIDlist[i]; forierr += SizeList[i] != Map.ElementSize(LIDlist[i]); } else { forierr += PIDlist[i] == MyPID; // If MyGID comes back false, the PID listed should be that of another proc } } EPETRA_TEST_ERR(forierr,ierr); delete [] MyGIDlist; delete [] GIDlist; delete [] PIDlist; delete [] LIDlist; delete [] SizeList; } delete [] MyGlobalElements1; delete [] MyElementSizeList; // Check RemoteIDList function (assumes all maps are linear, even if not stored that way) if (Map.LinearMap()) { int * GIDList = new int[3]; int * PIDList = new int[3]; int * LIDList = new int[3]; int MyPID = Map.Comm().MyPID(); int NumIDs = 0; //GIDList[NumIDs++] = Map.MaxAllGID()+1; // Should return -1 for both PID and LID if (Map.MinMyGID()-1>=Map.MinAllGID()) GIDList[NumIDs++] = Map.MinMyGID()-1; if (Map.MaxMyGID()+1<=Map.MaxAllGID()) GIDList[NumIDs++] = Map.MaxMyGID()+1; Map.RemoteIDList(NumIDs, GIDList, PIDList, LIDList); NumIDs = 0; //EPETRA_TEST_ERR(!(PIDList[NumIDs]==-1),ierr); //EPETRA_TEST_ERR(!(LIDList[NumIDs++]==-1),ierr); if (Map.MinMyGID()-1>=Map.MinAllGID()) EPETRA_TEST_ERR(!(PIDList[NumIDs++]==MyPID-1),ierr); if (Map.MaxMyGID()+1<=Map.MaxAllGID()) EPETRA_TEST_ERR(!(PIDList[NumIDs]==MyPID+1),ierr); if (Map.MaxMyGID()+1<=Map.MaxAllGID()) EPETRA_TEST_ERR(!(LIDList[NumIDs++]==0),ierr); delete [] GIDList; delete [] PIDList; delete [] LIDList; } return (ierr); }
void Epetra_Import::Construct( const Epetra_BlockMap & targetMap, const Epetra_BlockMap & sourceMap) { int i; // Build three ID lists: // NumSameIDs - Number of IDs in TargetMap and SourceMap that are identical, up to the first // nonidentical ID. // NumPermuteIDs - Number of IDs in SourceMap that must be indirectly loaded but are on this processor. // NumRemoteIDs - Number of IDs that are in SourceMap but not in TargetMap, and thus must be imported. int NumSourceIDs = sourceMap.NumMyElements(); int NumTargetIDs = targetMap.NumMyElements(); int_type *TargetGIDs = 0; if (NumTargetIDs>0) { TargetGIDs = new int_type[NumTargetIDs]; targetMap.MyGlobalElements(TargetGIDs); } int_type * SourceGIDs = 0; if (NumSourceIDs>0) { SourceGIDs = new int_type[NumSourceIDs]; sourceMap.MyGlobalElements(SourceGIDs); } int MinIDs = EPETRA_MIN(NumSourceIDs, NumTargetIDs); NumSameIDs_ = 0; for (i=0; i< MinIDs; i++) if (TargetGIDs[i]==SourceGIDs[i]) NumSameIDs_++; else break; // Find count of Target IDs that are truly remote and those that are local but permuted NumPermuteIDs_ = 0; NumRemoteIDs_ = 0; for (i=NumSameIDs_; i< NumTargetIDs; i++) if (sourceMap.MyGID(TargetGIDs[i])) NumPermuteIDs_++; // Check if Target GID is a local Source GID else NumRemoteIDs_++; // If not, then it is remote // Define remote and permutation lists int_type * RemoteGIDs=0; RemoteLIDs_ = 0; if (NumRemoteIDs_>0) { RemoteLIDs_ = new int[NumRemoteIDs_]; RemoteGIDs = new int_type[NumRemoteIDs_]; } if (NumPermuteIDs_>0) { PermuteToLIDs_ = new int[NumPermuteIDs_]; PermuteFromLIDs_ = new int[NumPermuteIDs_]; } NumPermuteIDs_ = 0; NumRemoteIDs_ = 0; for (i=NumSameIDs_; i< NumTargetIDs; i++) { if (sourceMap.MyGID(TargetGIDs[i])) { PermuteToLIDs_[NumPermuteIDs_] = i; PermuteFromLIDs_[NumPermuteIDs_++] = sourceMap.LID(TargetGIDs[i]); } else { //NumRecv_ +=TargetMap.ElementSize(i); // Count total number of entries to receive NumRecv_ +=targetMap.MaxElementSize(); // Count total number of entries to receive (currently need max) RemoteGIDs[NumRemoteIDs_] = TargetGIDs[i]; RemoteLIDs_[NumRemoteIDs_++] = i; } } if( NumRemoteIDs_>0 && !sourceMap.DistributedGlobal() ) ReportError("Warning in Epetra_Import: Serial Import has remote IDs. (Importing to Subset of Target Map)", 1); // Test for distributed cases int * RemotePIDs = 0; if (sourceMap.DistributedGlobal()) { if (NumRemoteIDs_>0) RemotePIDs = new int[NumRemoteIDs_]; int ierr = sourceMap.RemoteIDList(NumRemoteIDs_, RemoteGIDs, RemotePIDs, 0); // Get remote PIDs if (ierr) throw ReportError("Error in sourceMap.RemoteIDList call", ierr); //Get rid of IDs that don't exist in SourceMap if(NumRemoteIDs_>0) { int cnt = 0; for( i = 0; i < NumRemoteIDs_; ++i ) if( RemotePIDs[i] == -1 ) ++cnt; if( cnt ) { if( NumRemoteIDs_-cnt ) { int_type * NewRemoteGIDs = new int_type[NumRemoteIDs_-cnt]; int * NewRemotePIDs = new int[NumRemoteIDs_-cnt]; int * NewRemoteLIDs = new int[NumRemoteIDs_-cnt]; cnt = 0; for( i = 0; i < NumRemoteIDs_; ++i ) if( RemotePIDs[i] != -1 ) { NewRemoteGIDs[cnt] = RemoteGIDs[i]; NewRemotePIDs[cnt] = RemotePIDs[i]; NewRemoteLIDs[cnt] = targetMap.LID(RemoteGIDs[i]); ++cnt; } NumRemoteIDs_ = cnt; delete [] RemoteGIDs; delete [] RemotePIDs; delete [] RemoteLIDs_; RemoteGIDs = NewRemoteGIDs; RemotePIDs = NewRemotePIDs; RemoteLIDs_ = NewRemoteLIDs; ReportError("Warning in Epetra_Import: Target IDs not found in Source Map (Do you want to import to subset of Target Map?)", 1); } else { //valid RemoteIDs empty NumRemoteIDs_ = 0; delete [] RemoteGIDs; RemoteGIDs = 0; delete [] RemotePIDs; RemotePIDs = 0; } } } //Sort Remote IDs by processor so DoReverses will work Epetra_Util util; if(targetMap.GlobalIndicesLongLong()) { util.Sort(true,NumRemoteIDs_,RemotePIDs,0,0, 1,&RemoteLIDs_, 1,(long long**)&RemoteGIDs); } else if(targetMap.GlobalIndicesInt()) { int* ptrs[2] = {RemoteLIDs_, (int*)RemoteGIDs}; util.Sort(true,NumRemoteIDs_,RemotePIDs,0,0,2,&ptrs[0], 0, 0); } else { throw ReportError("Epetra_Import::Epetra_Import: GlobalIndices Internal Error", -1); } Distor_ = sourceMap.Comm().CreateDistributor(); // Construct list of exports that calling processor needs to send as a result // of everyone asking for what it needs to receive. bool Deterministic = true; int_type* tmp_ExportLIDs; //Export IDs come in as GIDs ierr = Distor_->CreateFromRecvs( NumRemoteIDs_, RemoteGIDs, RemotePIDs, Deterministic, NumExportIDs_, tmp_ExportLIDs, ExportPIDs_ ); if (ierr!=0) throw ReportError("Error in Epetra_Distributor.CreateFromRecvs()", ierr); // Export IDs come in as GIDs, convert to LIDs if(targetMap.GlobalIndicesLongLong()) { ExportLIDs_ = new int[NumExportIDs_]; for (i=0; i< NumExportIDs_; i++) { if (ExportPIDs_[i] < 0) throw ReportError("targetMap requested a GID that is not in the sourceMap.", -1); ExportLIDs_[i] = sourceMap.LID(tmp_ExportLIDs[i]); NumSend_ += sourceMap.MaxElementSize(); // Count total number of entries to send (currently need max) } delete[] tmp_ExportLIDs; } else if(targetMap.GlobalIndicesInt()) { for (i=0; i< NumExportIDs_; i++) { if (ExportPIDs_[i] < 0) throw ReportError("targetMap requested a GID that is not in the sourceMap.", -1); tmp_ExportLIDs[i] = sourceMap.LID(tmp_ExportLIDs[i]); NumSend_ += sourceMap.MaxElementSize(); // Count total number of entries to send (currently need max) } ExportLIDs_ = reinterpret_cast<int *>(tmp_ExportLIDs); // Can't reach here if tmp_ExportLIDs is long long. } else { throw ReportError("Epetra_Import::Epetra_Import: GlobalIndices Internal Error", -1); } } if( NumRemoteIDs_>0 ) delete [] RemoteGIDs; if( NumRemoteIDs_>0 ) delete [] RemotePIDs; if (NumTargetIDs>0) delete [] TargetGIDs; if (NumSourceIDs>0) delete [] SourceGIDs; return; }
void Epetra_Import::Construct_Expert( const Epetra_BlockMap & targetMap, const Epetra_BlockMap & sourceMap, int NumRemotePIDs, const int * UserRemotePIDs, const int & UserNumExportIDs, const int * UserExportLIDs, const int * UserExportPIDs) { int i,ierr; // Build three ID lists: // NumSameIDs - Number of IDs in TargetMap and SourceMap that are identical, up to the first // nonidentical ID. // NumPermuteIDs - Number of IDs in SourceMap that must be indirectly loaded but are on this processor. // NumRemoteIDs - Number of IDs that are in SourceMap but not in TargetMap, and thus must be imported. int NumSourceIDs = sourceMap.NumMyElements(); int NumTargetIDs = targetMap.NumMyElements(); int_type *TargetGIDs = 0; if (NumTargetIDs>0) { TargetGIDs = new int_type[NumTargetIDs]; targetMap.MyGlobalElements(TargetGIDs); } int_type * SourceGIDs = 0; if (NumSourceIDs>0) { SourceGIDs = new int_type[NumSourceIDs]; sourceMap.MyGlobalElements(SourceGIDs); } int MinIDs = EPETRA_MIN(NumSourceIDs, NumTargetIDs); NumSameIDs_ = 0; for (i=0; i< MinIDs; i++) if (TargetGIDs[i]==SourceGIDs[i]) NumSameIDs_++; else break; // Find count of Target IDs that are truly remote and those that are local but permuted NumPermuteIDs_ = 0; NumRemoteIDs_ = 0; for (i=NumSameIDs_; i< NumTargetIDs; i++) if (sourceMap.MyGID(TargetGIDs[i])) NumPermuteIDs_++; // Check if Target GID is a local Source GID else NumRemoteIDs_++; // If not, then it is remote // Define remote and permutation lists int_type * RemoteGIDs=0; RemoteLIDs_ = 0; if (NumRemoteIDs_>0) { RemoteLIDs_ = new int[NumRemoteIDs_]; RemoteGIDs = new int_type[NumRemoteIDs_]; } if (NumPermuteIDs_>0) { PermuteToLIDs_ = new int[NumPermuteIDs_]; PermuteFromLIDs_ = new int[NumPermuteIDs_]; } NumPermuteIDs_ = 0; NumRemoteIDs_ = 0; for (i=NumSameIDs_; i< NumTargetIDs; i++) { if (sourceMap.MyGID(TargetGIDs[i])) { PermuteToLIDs_[NumPermuteIDs_] = i; PermuteFromLIDs_[NumPermuteIDs_++] = sourceMap.LID(TargetGIDs[i]); } else { //NumRecv_ +=TargetMap.ElementSize(i); // Count total number of entries to receive NumRecv_ +=targetMap.MaxElementSize(); // Count total number of entries to receive (currently need max) RemoteGIDs[NumRemoteIDs_] = TargetGIDs[i]; RemoteLIDs_[NumRemoteIDs_++] = i; } } if( NumRemoteIDs_>0 && !sourceMap.DistributedGlobal() ) ReportError("Warning in Epetra_Import: Serial Import has remote IDs. (Importing to Subset of Target Map)", 1); // Test for distributed cases int * RemotePIDs = 0; if (sourceMap.DistributedGlobal()) { if (NumRemoteIDs_>0) RemotePIDs = new int[NumRemoteIDs_]; #ifdef EPETRA_ENABLE_DEBUG int myeq = (NumRemotePIDs==NumRemoteIDs_); int globaleq=0; sourceMap.Comm().MinAll(&myeq,&globaleq,1); if(globaleq!=1) { printf("[%d] UserRemotePIDs count wrong %d != %d\n",sourceMap.Comm().MyPID(),NumRemotePIDs,NumRemoteIDs_); fflush(stdout); sourceMap.Comm().Barrier(); sourceMap.Comm().Barrier(); sourceMap.Comm().Barrier(); throw ReportError("Epetra_Import: UserRemotePIDs count wrong"); } #endif if(NumRemotePIDs==NumRemoteIDs_){ // Since I need to sort these, I'll copy them for(i=0; i<NumRemoteIDs_; i++) RemotePIDs[i] = UserRemotePIDs[i]; } //Get rid of IDs that don't exist in SourceMap if(NumRemoteIDs_>0) { int cnt = 0; for( i = 0; i < NumRemoteIDs_; ++i ) if( RemotePIDs[i] == -1 ) ++cnt; if( cnt ) { if( NumRemoteIDs_-cnt ) { int_type * NewRemoteGIDs = new int_type[NumRemoteIDs_-cnt]; int * NewRemotePIDs = new int[NumRemoteIDs_-cnt]; int * NewRemoteLIDs = new int[NumRemoteIDs_-cnt]; cnt = 0; for( i = 0; i < NumRemoteIDs_; ++i ) if( RemotePIDs[i] != -1 ) { NewRemoteGIDs[cnt] = RemoteGIDs[i]; NewRemotePIDs[cnt] = RemotePIDs[i]; NewRemoteLIDs[cnt] = targetMap.LID(RemoteGIDs[i]); ++cnt; } NumRemoteIDs_ = cnt; delete [] RemoteGIDs; delete [] RemotePIDs; delete [] RemoteLIDs_; RemoteGIDs = NewRemoteGIDs; RemotePIDs = NewRemotePIDs; RemoteLIDs_ = NewRemoteLIDs; ReportError("Warning in Epetra_Import: Target IDs not found in Source Map (Do you want to import to subset of Target Map?)", 1); } else { //valid RemoteIDs empty NumRemoteIDs_ = 0; delete [] RemoteGIDs; RemoteGIDs = 0; delete [] RemotePIDs; RemotePIDs = 0; } } } //Sort Remote IDs by processor so DoReverses will work Epetra_Util util; if(targetMap.GlobalIndicesLongLong()) { util.Sort(true,NumRemoteIDs_,RemotePIDs,0,0, 1,&RemoteLIDs_, 1,(long long**)&RemoteGIDs); } else if(targetMap.GlobalIndicesInt()) { int* ptrs[2] = {RemoteLIDs_, (int*)RemoteGIDs}; util.Sort(true,NumRemoteIDs_,RemotePIDs,0,0,2,&ptrs[0], 0, 0); } else { throw ReportError("Epetra_Import::Epetra_Import: GlobalIndices Internal Error", -1); } // Build distributor & Export lists Distor_ = sourceMap.Comm().CreateDistributor(); NumExportIDs_=UserNumExportIDs; ExportLIDs_ = new int[NumExportIDs_]; ExportPIDs_ = new int[NumExportIDs_]; for(i=0; i<NumExportIDs_; i++) { ExportPIDs_[i] = UserExportPIDs[i]; ExportLIDs_[i] = UserExportLIDs[i]; } #ifdef HAVE_MPI Epetra_MpiDistributor* MpiDistor = dynamic_cast< Epetra_MpiDistributor*>(Distor_); bool Deterministic = true; if(MpiDistor) ierr=MpiDistor->CreateFromSendsAndRecvs(NumExportIDs_,ExportPIDs_, NumRemoteIDs_, RemoteGIDs, RemotePIDs,Deterministic); else ierr=-10; #else ierr=-20; #endif if (ierr!=0) throw ReportError("Error in Epetra_Distributor.CreateFromRecvs()", ierr); } if( NumRemoteIDs_>0 ) delete [] RemoteGIDs; if( NumRemoteIDs_>0 ) delete [] RemotePIDs; if (NumTargetIDs>0) delete [] TargetGIDs; if (NumSourceIDs>0) delete [] SourceGIDs; #ifdef EPETRA_ENABLE_DEBUG // Sanity check to make sure we got the import right Epetra_IntVector Source(sourceMap); Epetra_IntVector Target(targetMap); for(i=0; i<Source.MyLength(); i++) Source[i] = (int) (Source.Map().GID(i) % INT_MAX); Target.PutValue(-1); Target.Import(Source,*this,Insert); bool test_passed=true; for(i=0; i<Target.MyLength(); i++){ if(Target[i] != Target.Map().GID(i) % INT_MAX) test_passed=false; } if(!test_passed) { printf("[%d] PROCESSOR has a mismatch... prepearing to crash or hang!\n",sourceMap.Comm().MyPID()); fflush(stdout); sourceMap.Comm().Barrier(); sourceMap.Comm().Barrier(); sourceMap.Comm().Barrier(); throw ReportError("Epetra_Import: ERROR. User provided IDs do not match what an import generates."); } #endif return; }