54 #include "../epetra_test_err.h"
64 int main(
int argc,
char *argv[]) {
68 if ((argv[1][0] ==
'-') && (argv[1][1] ==
'v')) {
78 MPI_Init(&argc,&argv);
87 int MyPID = Comm.
MyPID();
89 int verbose_int = verbose ? 1 : 0;
91 verbose = verbose_int==1 ?
true :
false;
93 if (verbose && MyPID==0)
111 if (returnierr == 0) {
112 cout <<
"Epetra_Directory tests passed."<<endl;
115 cout <<
"Epetra_Directory tests failed."<<endl;
127 int myPID = Comm.
MyPID();
130 if (numProcs < 2)
return(0);
132 int myFirstID = (myPID+1)*(myPID+1);
133 int myNumIDs = 3+myPID;
135 long long* myIDs =
new long long[myNumIDs];
137 for(i=0; i<myNumIDs; ++i) {
138 myIDs[i] = myFirstID+i;
146 if (proc >= numProcs) proc = 0;
148 int procNumIDs = 3+proc;
149 long long procFirstID = (
long long)(proc+1)*(proc+1);
150 long long procLastID = procFirstID+procNumIDs - 1;
156 blkmap, 1, &procFirstID,
157 &queryProc1, NULL, NULL
160 blkmap, 1, &procLastID,
161 &queryProc2, NULL, NULL
166 if (queryProc1 != proc || queryProc2 != proc) {
181 int myPID = Comm.
MyPID();
184 if (numProcs < 2)
return(0);
186 int myFirstID = (numProcs-myPID)*(numProcs-myPID);
189 long long* myIDs =
new long long[myNumIDs];
191 for(i=0; i<myNumIDs; ++i) {
192 myIDs[i] = myFirstID+i;
200 if (proc >= numProcs) proc = 0;
203 long long procFirstID = (
long long)(numProcs-proc)*(numProcs-proc);
204 long long procLastID = procFirstID+procNumIDs - 1;
210 &queryProc1, NULL, NULL);
212 &queryProc2, NULL, NULL);
216 if (queryProc1 != proc || queryProc2 != proc) {
228 int myPID = Comm.
MyPID();
231 if (numProcs < 2)
return(0);
233 int myFirstID = (myPID+1)*(myPID+1);
236 long long* myIDs =
new long long[myNumIDs];
238 for(i=0; i<myNumIDs-1; ++i) {
239 myIDs[i] = myFirstID+i;
242 int nextProc = myPID+1;
243 if (nextProc >= numProcs) nextProc = 0;
245 int nextProcFirstID = (nextProc+1)*(nextProc+1);
246 myIDs[myNumIDs-1] = nextProcFirstID;
266 int myPID = Comm.
MyPID();
269 if (numProcs < 2)
return(0);
273 int numMyGIDs = 2*num;
274 int myFirstGID = myPID*num;
276 long long* myGIDs =
new long long[numMyGIDs];
278 for(
int i=0; i<numMyGIDs; ++i) {
279 myGIDs[i] = myFirstGID+i;
282 Epetra_Map overlappingmap((
long long)-1, numMyGIDs, myGIDs, 0, Comm);
286 long long numGlobal0 = overlappingmap.NumGlobalElements64();
288 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES // FIXME
292 bool use_high_sharing_proc =
true;
301 if (numGlobal1 != numGlobal2) {
307 if (numGlobal0 <= numGlobal1) {
316 if ((myPID==0 || myPID==numProcs-1) && numLocal1 == numLocal2) {
326 int myPID = Comm.
MyPID();
329 if (numProcs < 2)
return(0);
333 int numMyGIDs = 2*num;
334 int myFirstGID = myPID*num;
336 long long* myGIDs =
new long long[numMyGIDs];
337 int* sizes =
new int[numMyGIDs];
339 for(
int i=0; i<numMyGIDs; ++i) {
340 myGIDs[i] = myFirstGID+i;
341 sizes[i] = myFirstGID+i+1;
344 Epetra_BlockMap overlappingmap((
long long)-1, numMyGIDs, myGIDs, sizes, 0, Comm);
349 long long numGlobal0 = overlappingmap.NumGlobalElements64();
351 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES // FIXME
355 bool use_high_sharing_proc =
true;
364 if (numGlobal1 != numGlobal2) {
370 if (numGlobal0 <= numGlobal1) {
379 if ((myPID==0 || myPID==numProcs-1) && numLocal1 == numLocal2) {
Epetra_Map: A class for partitioning vectors and matrices.
int directory_test_3(Epetra_Comm &Comm)
static Epetra_Map Create_OneToOne_Map(const Epetra_Map &usermap, bool high_rank_proc_owns_shared=false)
Epetra_Util Create_OneToOne_Map function.
#define EPETRA_TEST_ERR(a, b)
long long NumGlobalElements64() const
virtual Epetra_Directory * CreateDirectory(const Epetra_BlockMap &Map) const =0
Create a directory object for the given Epetra_BlockMap.
static void SetTracebackMode(int TracebackModeValue)
Set the value of the Epetra_Object error traceback report mode.
virtual bool GIDsAllUniquelyOwned() const =0
GIDsAllUniquelyOwned: returns true if all GIDs appear on just one processor.
int MyPID() const
Return my process ID.
Epetra_MpiComm: The Epetra MPI Communication Class.
virtual int GetDirectoryEntries(const Epetra_BlockMap &Map, const int NumEntries, const int *GlobalEntries, int *Procs, int *LocalEntries, int *EntrySizes, bool high_rank_sharing_procs=false) const =0
GetDirectoryEntries : Returns proc and local id info for non-local map entries.
std::string Epetra_Version()
int directory_test_4(Epetra_Comm &Comm)
virtual int MyPID() const =0
Return my process ID.
int NumMyElements() const
Number of elements on the calling processor.
Epetra_Directory: This class is a pure virtual class whose interface allows Epetra_Map and Epetr_Bloc...
static Epetra_BlockMap Create_OneToOne_BlockMap(const Epetra_BlockMap &usermap, bool high_rank_proc_owns_shared=false)
Epetra_Util Create_OneToOne_Map function.
Epetra_Comm: The Epetra Communication Abstract Base Class.
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
int directory_test_5(Epetra_Comm &Comm)
int directory_test_2(Epetra_Comm &Comm)
Epetra_SerialComm: The Epetra Serial Communication Class.
int directory_test_1(Epetra_Comm &Comm)
virtual int NumProc() const =0
Returns total number of processes.
int main(int argc, char *argv[])
int Broadcast(double *MyVals, int Count, int Root) const
Epetra_SerialComm Broadcast function.