Epetra Package Browser (Single Doxygen Collection) Development
Loading...
Searching...
No Matches
Epetra_BlockMap.h
Go to the documentation of this file.
1/*
2//@HEADER
3// ************************************************************************
4//
5// Epetra: Linear Algebra Services Package
6// Copyright 2011 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39//
40// ************************************************************************
41//@HEADER
42*/
43
44#ifndef EPETRA_BLOCKMAP_H
45#define EPETRA_BLOCKMAP_H
46
47#include "Epetra_ConfigDefs.h"
48#include "Epetra_Object.h"
49#include "Epetra_BlockMapData.h"
50
51
53
194class EPETRA_LIB_DLL_EXPORT Epetra_BlockMap: public Epetra_Object {
195 friend class Epetra_Directory;
196 friend class Epetra_LocalMap;
197 public:
199
200
224#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
225 Epetra_BlockMap(int NumGlobalElements, int ElementSize, int IndexBase, const Epetra_Comm& Comm);
226#endif
227#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
228 Epetra_BlockMap(long long NumGlobalElements, int ElementSize, int IndexBase, const Epetra_Comm& Comm);
229 Epetra_BlockMap(long long NumGlobalElements, int ElementSize, long long IndexBase, const Epetra_Comm& Comm);
230#endif
231
233
262#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
263 Epetra_BlockMap(int NumGlobalElements, int NumMyElements,
264 int ElementSize, int IndexBase, const Epetra_Comm& Comm);
265#endif
266#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
267 Epetra_BlockMap(long long NumGlobalElements, int NumMyElements,
268 int ElementSize, int IndexBase, const Epetra_Comm& Comm);
269 Epetra_BlockMap(long long NumGlobalElements, int NumMyElements,
270 int ElementSize, long long IndexBase, const Epetra_Comm& Comm);
271#endif
272
274
310#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
311 Epetra_BlockMap(int NumGlobalElements, int NumMyElements,
312 const int *MyGlobalElements,
313 int ElementSize, int IndexBase, const Epetra_Comm& Comm);
314#endif
315#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
316 Epetra_BlockMap(long long NumGlobalElements, int NumMyElements,
317 const long long *MyGlobalElements,
318 int ElementSize, int IndexBase, const Epetra_Comm& Comm);
319 Epetra_BlockMap(long long NumGlobalElements, int NumMyElements,
320 const long long *MyGlobalElements,
321 int ElementSize, long long IndexBase, const Epetra_Comm& Comm);
322#endif
323
325
362#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
363 Epetra_BlockMap(int NumGlobalElements, int NumMyElements,
364 const int *MyGlobalElements,
365 const int *ElementSizeList, int IndexBase,
366 const Epetra_Comm& Comm);
367#endif
368#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
369 Epetra_BlockMap(long long NumGlobalElements, int NumMyElements,
370 const long long *MyGlobalElements,
371 const int *ElementSizeList, int IndexBase,
372 const Epetra_Comm& Comm);
373 Epetra_BlockMap(long long NumGlobalElements, int NumMyElements,
374 const long long *MyGlobalElements,
375 const int *ElementSizeList, long long IndexBase,
376 const Epetra_Comm& Comm);
377#endif
378
379#if defined(EPETRA_NO_32BIT_GLOBAL_INDICES) && defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
380 // default implementation so that no compiler/linker error in case neither 32 nor 64
381 // bit indices present.
382 Epetra_BlockMap() {}
383#endif
385
387#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
388 Epetra_BlockMap(long long NumGlobal_Elements, int NumMy_Elements,
389 const long long * myGlobalElements,
390 int ElementSize, int indexBase,
391 const Epetra_Comm& comm,
392 bool UserIsDistributedGlobal,
393 long long UserMinAllGID, long long UserMaxAllGID);
394 Epetra_BlockMap(long long NumGlobal_Elements, int NumMy_Elements,
395 const long long * myGlobalElements,
396 int ElementSize, long long indexBase,
397 const Epetra_Comm& comm,
398 bool UserIsDistributedGlobal,
399 long long UserMinAllGID, long long UserMaxAllGID);
400#endif
401
402#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
403 Epetra_BlockMap(int NumGlobal_Elements, int NumMy_Elements,
404 const int * myGlobalElements,
405 int ElementSize, int indexBase,
406 const Epetra_Comm& comm,
407 bool UserIsDistributedGlobal,
408 int UserMinAllGID, int UserMaxAllGID);
409#endif
410
411
414
416 virtual ~Epetra_BlockMap(void);
418
420
421
429#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
430 int RemoteIDList(int NumIDs, const int * GIDList, int * PIDList, int * LIDList) const {
431 return(RemoteIDList(NumIDs, GIDList, PIDList, LIDList, 0));
432 };
433#endif
434#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
435 int RemoteIDList(int NumIDs, const long long * GIDList, int * PIDList, int * LIDList) const {
436 return(RemoteIDList(NumIDs, GIDList, PIDList, LIDList, 0));
437 };
438#endif
439
441
446#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
447 int RemoteIDList(int NumIDs, const int * GIDList, int * PIDList, int * LIDList, int * SizeList) const;
448#endif
449#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
450 int RemoteIDList(int NumIDs, const long long * GIDList, int * PIDList, int * LIDList, int * SizeList) const;
451#endif
452
454#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
455 int LID(int GID) const;
456#endif
457#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
458 int LID(long long GID) const;
459#endif
460
461#if defined(EPETRA_NO_32BIT_GLOBAL_INDICES) && defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
462 // default implementation so that no compiler/linker error in case neither 32 nor 64
463 // bit indices present.
464 int LID(long long GID) const { return -1; }
465 bool MyGID(long long GID_in) const { return false; }
466#endif
467
469#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
470 int GID(int LID) const;
471#endif
472 long long GID64(int LID) const;
473
475 int FindLocalElementID(int PointID, int & ElementID, int & ElementOffset) const;
476
478#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
479 bool MyGID(int GID_in) const {return(LID(GID_in)!=-1);};
480#endif
481#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
482 bool MyGID(long long GID_in) const {return(LID(GID_in)!=-1);};
483#endif
484
486 // bool MyLID(int LID_in) const {return(GID64(LID_in)!=BlockMapData_->IndexBase_-1);};
487 bool MyLID(int lid) const {
488 if ((BlockMapData_->NumMyElements_ == 0) ||
489 (lid < BlockMapData_->MinLID_) || (lid > BlockMapData_->MaxLID_)) {
490 return false;
491 }
492 return true;
493 }
494
496#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
497 int MinAllGID() const {
498 if(GlobalIndicesInt())
499 return (int) MinAllGID64();
500 throw "Epetra_BlockMap::MinAllGID: GlobalIndices not int.";
501 }
502#endif
503 long long MinAllGID64() const {return(BlockMapData_->MinAllGID_);};
504
506#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
507 int MaxAllGID() const {
508 if(GlobalIndicesInt())
509 return (int) MaxAllGID64();
510 throw "Epetra_BlockMap::MaxAllGID: GlobalIndices not int.";
511 }
512#endif
513 long long MaxAllGID64() const {return(BlockMapData_->MaxAllGID_);};
514
516#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
517 int MinMyGID() const {
518 if(GlobalIndicesInt())
519 return (int) MinMyGID64();
520 throw "Epetra_BlockMap::MinMyGID: GlobalIndices not int.";
521 }
522#endif
523 long long MinMyGID64() const {return(BlockMapData_->MinMyGID_);};
524
526#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
527 int MaxMyGID() const {
528 if(GlobalIndicesInt())
529 return (int) MaxMyGID64();
530 throw "Epetra_BlockMap::MaxMyGID: GlobalIndices not int.";
531 }
532#endif
533 long long MaxMyGID64() const {return(BlockMapData_->MaxMyGID_);};
534
536 int MinLID() const {return(BlockMapData_->MinLID_);};
537
539 int MaxLID() const {return(BlockMapData_->MaxLID_);};
541
543
544
545#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
546 int NumGlobalElements() const {
547 if(GlobalIndicesInt())
548 return (int) NumGlobalElements64();
549 throw "Epetra_BlockMap::NumGlobalElements: GlobalIndices not int.";
550 }
551#endif
552 long long NumGlobalElements64() const {return(BlockMapData_->NumGlobalElements_);};
553
555 int NumMyElements() const {return(BlockMapData_->NumMyElements_);};
556
558#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
559 int MyGlobalElements(int * MyGlobalElementList) const;
560#endif
561#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
562 int MyGlobalElements(long long * MyGlobalElementList) const;
563#endif
564
565#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
566 int MyGlobalElementsPtr(int *& MyGlobalElementList) const;
567#endif
568#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
569 int MyGlobalElementsPtr(long long *& MyGlobalElementList) const;
570#endif
571
573 int ElementSize() const {return(BlockMapData_->ElementSize_);};
574
576 int ElementSize(int LID) const;
577
579
582 int FirstPointInElement(int LID) const;
583
584#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
586 int IndexBase() const {
587 if(GlobalIndicesInt() || IndexBase64() == (long long) static_cast<int>(IndexBase64()))
588 return (int) IndexBase64();
589 throw "Epetra_BlockMap::IndexBase: GlobalIndices not int and IndexBase cannot fit an int.";
590 }
591#endif
592 long long IndexBase64() const {return(BlockMapData_->IndexBase_);};
593
595#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
596 int NumGlobalPoints() const {
597 if(GlobalIndicesInt())
598 return (int) NumGlobalPoints64();
599 throw "Epetra_BlockMap::NumGlobalPoints: GlobalIndices not int.";
600 }
601#endif
602 long long NumGlobalPoints64() const {return(BlockMapData_->NumGlobalPoints_);};
603
605 int NumMyPoints() const {return(BlockMapData_->NumMyPoints_);};
606
608 int MinMyElementSize() const {return(BlockMapData_->MinMyElementSize_);};
609
611 int MaxMyElementSize() const {return(BlockMapData_->MaxMyElementSize_);};
612
614 int MinElementSize() const {return(BlockMapData_->MinElementSize_);};
615
617 int MaxElementSize() const {return(BlockMapData_->MaxElementSize_);};
619
621
622
628 bool UniqueGIDs() const {return(IsOneToOne());};
629
630/*
631*******************************************************************************
632 Ideally GlobalIndicesInt and GlobalIndicesLongLong should be within the
633 preprocessor macros and any code using them should also be within the
634 corresponding macro. However, when initially moving to 64-bit we did not
635 have macros and all the code is written using run-time checks. In future,
636 the code can be converted to follow the macro system. Hence this comment.
637 -- Chetan Jhurani
638
639 Future code:
640
641#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
643 bool GlobalIndicesInt() const { return BlockMapData_->GlobalIndicesInt_; }
644#endif
645#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
647 bool GlobalIndicesLongLong() const { return BlockMapData_->GlobalIndicesLongLong_; }
648#endif
649*******************************************************************************
650*/
651
653 bool GlobalIndicesInt() const { return BlockMapData_->GlobalIndicesInt_; }
655 bool GlobalIndicesLongLong() const { return BlockMapData_->GlobalIndicesLongLong_; }
656
657 template<typename int_type>
659
660 bool GlobalIndicesTypeValid() const { return BlockMapData_->GlobalIndicesInt_ || BlockMapData_->GlobalIndicesLongLong_; }
661
663 {
664 return
665#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
666 GlobalIndicesInt() == other.GlobalIndicesInt() &&
667#endif
668#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
669 GlobalIndicesLongLong() == other.GlobalIndicesLongLong() &&
670#endif
671 true;
672 }
673
675 bool ConstantElementSize() const {return(BlockMapData_->ConstantElementSize_);};
676
679 bool SameBlockMapDataAs(const Epetra_BlockMap & Map) const;
680
682 bool SameAs(const Epetra_BlockMap & Map) const;
683
685
688 bool PointSameAs(const Epetra_BlockMap & Map) const;
689
691 bool LinearMap() const {return(BlockMapData_->LinearMap_);};
692
694 bool DistributedGlobal() const {return(BlockMapData_->DistributedGlobal_);};
696
698
699
701#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
702 int * MyGlobalElements() const;
703#endif
704#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
705 long long * MyGlobalElements64() const;
706#endif
707
708 // Helper function to avoid scattering ifdef in other code.
709 void MyGlobalElements(const int*& IntGIDs, const long long*& LLGIDs) const {
710#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
711 if(GlobalIndicesInt()) {
712 IntGIDs = MyGlobalElements();
713 return;
714 }
715#endif
716#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
717 if(GlobalIndicesLongLong()) {
718 LLGIDs = MyGlobalElements64();
719 return;
720 }
721#endif
722 }
723
724 // Helper function to avoid scattering ifdef in other code.
725 void MyGlobalElements(int*& IntGIDs, long long*& LLGIDs) {
726#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
727 if(GlobalIndicesInt()) {
728 IntGIDs = MyGlobalElements();
729 return;
730 }
731#endif
732#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
733 if(GlobalIndicesLongLong()) {
734 LLGIDs = MyGlobalElements64();
735 return;
736 }
737#endif
738 }
739
741
744 int * FirstPointInElementList() const;
745
747 int * ElementSizeList() const;
748
750 int * PointToElementList() const;
751
753 int ElementSizeList(int * ElementSizeList)const;
754
756 int FirstPointInElementList(int * FirstPointInElementList)const;
757
759 int PointToElementList(int * PointToElementList) const;
760
762
764
765
767 virtual void Print(std::ostream & os) const;
768
770 const Epetra_Comm & Comm() const {return(*BlockMapData_->Comm_);}
771
772 bool IsOneToOne() const;
773
776
778
780
781
783
784 int ReferenceCount() const {return(BlockMapData_->ReferenceCount());}
785
787
788 const Epetra_BlockMapData * DataPtr() const {return(BlockMapData_);}
789
833 Epetra_BlockMap * RemoveEmptyProcesses() const;
834
862 Epetra_BlockMap* ReplaceCommWithSubset(const Epetra_Comm * Comm) const;
863
865
866 private: // These need to be accessible to derived map classes.
867
868 void GlobalToLocalSetup();
869 bool DetermineIsOneToOne() const;
870 bool IsDistributedGlobal(long long NumGlobalElements, int NumMyElements) const;
871 void CheckValidNGE(long long NumGlobalElements);
872 void EndOfConstructorOps();
873
874 protected:
875 void CleanupData();
876
878
879private:
880
881
882 void ConstructAutoUniform(long long NumGlobal_Elements, int Element_Size,
883 long long Index_Base, const Epetra_Comm& comm, bool IsLongLong);
884
885 void ConstructUserLinear(long long NumGlobal_Elements, int NumMy_Elements,
886 int Element_Size, long long Index_Base, const Epetra_Comm& comm, bool IsLongLong);
887
888 template<typename int_type>
889 void ConstructUserConstant(int_type NumGlobal_Elements, int NumMy_Elements,
890 const int_type * myGlobalElements,
891 int Element_Size, int_type indexBase,
892 const Epetra_Comm& comm, bool IsLongLong);
893
894 template<typename int_type>
895 void ConstructUserVariable(int_type NumGlobal_Elements, int NumMy_Elements,
896 const int_type * myGlobalElements,
897 const int *elementSizeList, int_type indexBase,
898 const Epetra_Comm& comm, bool IsLongLong);
899
900 template<typename int_type>
901 void ConstructUserConstantNoComm(int_type NumGlobal_Elements, int NumMy_Elements,
902 const int_type * myGlobalElements,
903 int ElementSize, int_type indexBase,
904 const Epetra_Comm& comm, bool IsLongLong,
905 bool UserIsDistributedGlobal,
906 int_type UserMinAllGID, int_type UserMaxAllGID);
907
908 template<typename int_type>
909 int_type& MyGlobalElementVal(int i);
910
911 template<typename int_type>
912 int_type MyGlobalElementValGet(int i);
913
914 template<typename int_type>
916
917 template<typename int_type>
918 void TGlobalToLocalSetup();
919};
920
921#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
922template<> inline bool Epetra_BlockMap::GlobalIndicesIsType<long long>() const { return BlockMapData_->GlobalIndicesLongLong_; }
923template<> inline long long& Epetra_BlockMap::MyGlobalElementVal<long long>(int i) { return BlockMapData_->MyGlobalElements_LL_[i]; }
924template<> inline long long Epetra_BlockMap::MyGlobalElementValGet<long long>(int i) { return BlockMapData_->MyGlobalElements_LL_[i]; }
925template<> inline int Epetra_BlockMap::SizeMyGlobalElement<long long>(int n) { return BlockMapData_->MyGlobalElements_LL_.Size(n); }
926#endif
927
928#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
929template<> inline bool Epetra_BlockMap::GlobalIndicesIsType<int>() const { return BlockMapData_->GlobalIndicesInt_; }
930template<> inline int& Epetra_BlockMap::MyGlobalElementVal<int> (int i) { return BlockMapData_->MyGlobalElements_int_[i]; }
931template<> inline int Epetra_BlockMap::MyGlobalElementValGet<int> (int i) { return BlockMapData_->MyGlobalElements_int_[i]; }
932template<> inline int Epetra_BlockMap::SizeMyGlobalElement<int> (int n) { return BlockMapData_->MyGlobalElements_int_.Size(n); }
933#endif
934
935#endif /* EPETRA_BLOCKMAP_H */
Epetra_BlockMapData: The Epetra BlockMap Data Class.
Epetra_IntSerialDenseVector MyGlobalElements_int_
Epetra_LongLongSerialDenseVector MyGlobalElements_LL_
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
long long IndexBase64() const
int MinAllGID() const
Returns the minimum global ID across the entire map.
int SizeMyGlobalElement(int n)
int_type & MyGlobalElementVal(int i)
int RemoteIDList(int NumIDs, const int *GIDList, int *PIDList, int *LIDList) const
Returns the processor IDs and corresponding local index value for a given list of global indices.
int MaxElementSize() const
Maximum element size across all processors.
bool GlobalIndicesIsType() const
int MinMyElementSize() const
Minimum element size on the calling processor.
int_type MyGlobalElementValGet(int i)
long long MaxAllGID64() const
int MinMyGID() const
Returns the minimum global ID owned by this processor.
long long MaxMyGID64() const
bool GlobalIndicesTypeMatch(const Epetra_BlockMap &other) const
bool GlobalIndicesTypeValid() const
bool MyGID(long long GID_in) const
bool DistributedGlobal() const
Returns true if map is defined across more than one processor.
int NumMyPoints() const
Number of local points for this map; equals the sum of all element sizes on the calling processor.
int ReferenceCount() const
Returns the reference count of BlockMapData.
bool LinearMap() const
Returns true if the global ID space is contiguously divided (but not necessarily uniformly) across al...
int IndexBase() const
Index base for this map.
Epetra_BlockMapData * BlockMapData_
int RemoteIDList(int NumIDs, const long long *GIDList, int *PIDList, int *LIDList) const
int MaxAllGID() const
Returns the maximum global ID across the entire map.
int ElementSize() const
Returns the size of elements in the map; only valid if map has constant element size.
void MyGlobalElements(int *&IntGIDs, long long *&LLGIDs)
int MaxMyElementSize() const
Maximum element size on the calling processor.
int MinLID() const
The minimum local index value on the calling processor.
bool UniqueGIDs() const
Returns true if map GIDs are 1-to-1.
long long NumGlobalElements64() const
void MyGlobalElements(const int *&IntGIDs, const long long *&LLGIDs) const
bool GlobalIndicesInt() const
Returns true if map create with int NumGlobalElements.
int MinElementSize() const
Minimum element size across all processors.
int NumGlobalElements() const
Number of elements across all processors.
const Epetra_Comm & Comm() const
Access function for Epetra_Comm communicator.
bool ConstantElementSize() const
Returns true if map has constant element size.
const Epetra_BlockMapData * DataPtr() const
Returns a pointer to the BlockMapData instance this BlockMap uses.
long long MinAllGID64() const
long long MinMyGID64() const
int NumMyElements() const
Number of elements on the calling processor.
bool MyLID(int lid) const
Returns true if the LID passed in belongs to the calling processor in this map, otherwise returns fal...
long long NumGlobalPoints64() const
int MaxLID() const
The maximum local index value on the calling processor.
int NumGlobalPoints() const
Number of global points for this map; equals the sum of all element sizes across all processors.
bool MyGID(int GID_in) const
Returns true if the GID passed in belongs to the calling processor in this map, otherwise returns fal...
bool GlobalIndicesLongLong() const
Returns true if map create with long long NumGlobalElements.
int MaxMyGID() const
Returns the maximum global ID owned by this processor.
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:73
Epetra_Directory: This class is a pure virtual class whose interface allows Epetra_Map and Epetr_Bloc...
int Size(int Length_in)
Set length of a Epetra_IntSerialDenseVector object; init values to zero.
Epetra_LocalMap: A class for replicating vectors and matrices across multiple processors.
int Size(int Length_in)
Set length of a Epetra_LongLongSerialDenseVector object; init values to zero.
Epetra_Object: The base Epetra class.
Definition: Epetra_Object.h:57
virtual void Print(std::ostream &os) const
Print object to an output stream Print method.
Epetra_Object & operator=(const Epetra_Object &src)