Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_IntMultiVector.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_INTMULTIVECTOR_H
45 #define EPETRA_INTMULTIVECTOR_H
46 
47 #if defined(Epetra_SHOW_DEPRECATED_WARNINGS)
48 #ifdef __GNUC__
49 #warning "The Epetra package is deprecated"
50 #endif
51 #endif
52 
53 
54 
55 class Epetra_Comm;
56 class Epetra_BlockMap;
57 class Epetra_Map;
58 class Epetra_Import;
59 class Epetra_Export;
60 class Epetra_Distributor;
61 class Epetra_IntVector;
62 
63 #include "Epetra_ConfigDefs.h"
64 #include "Epetra_DistObject.h"
65 #include "Epetra_CompObject.h"
66 #include "Epetra_BLAS.h"
67 #include "Epetra_Util.h"
68 
70 
186 //==========================================================================
187 class EPETRA_LIB_DLL_EXPORT Epetra_IntMultiVector: public Epetra_DistObject, public Epetra_CompObject, public Epetra_BLAS {
188 
189  public:
190 
192 
193 
211  Epetra_IntMultiVector(const Epetra_BlockMap& Map, int NumVectors, bool zeroOut = true);
212 
214 
216 
218 
238  int *A, int MyLDA, int NumVectors);
239 
241 
257  int **ArrayOfPointers, int NumVectors);
258 
260 
275  const Epetra_IntMultiVector& Source, int *Indices, int NumVectors);
276 
278 
293  const Epetra_IntMultiVector& Source, int StartIndex,
294  int NumVectors);
295 
297  virtual ~Epetra_IntMultiVector();
299 
301 
302 
304 
323 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
324  int ReplaceGlobalValue(int GlobalRow, int VectorIndex, int OrdinalValue);
325 #endif
326 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
327  int ReplaceGlobalValue(long long GlobalRow, int VectorIndex, int OrdinalValue);
328 #endif
329 
330 
332 
350 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
351  int ReplaceGlobalValue(int GlobalBlockRow, int BlockRowOffset, int VectorIndex, int OrdinalValue);
352 #endif
353 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
354  int ReplaceGlobalValue(long long GlobalBlockRow, int BlockRowOffset, int VectorIndex, int OrdinalValue);
355 #endif
356 
357 
359 
378 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
379  int SumIntoGlobalValue(int GlobalRow, int VectorIndex, int OrdinalValue);
380 #endif
381 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
382  int SumIntoGlobalValue(long long GlobalRow, int VectorIndex, int OrdinalValue);
383 #endif
384 
385 
387 
405 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
406  int SumIntoGlobalValue(int GlobalBlockRow, int BlockRowOffset, int VectorIndex, int OrdinalValue);
407 #endif
408 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
409  int SumIntoGlobalValue(long long GlobalBlockRow, int BlockRowOffset, int VectorIndex, int OrdinalValue);
410 #endif
411 
413 
437  int ReplaceMyValue(int MyRow, int VectorIndex, int OrdinalValue);
438 
439 
441 
459  int ReplaceMyValue(int MyBlockRow, int BlockRowOffset, int VectorIndex, int OrdinalValue);
460 
461 
463 
482  int SumIntoMyValue(int MyRow, int VectorIndex, int OrdinalValue);
483 
484 
486 
504  int SumIntoMyValue(int MyBlockRow, int BlockRowOffset, int VectorIndex, int OrdinalValue);
505 
507 
513  int PutScalar (int OrdinalConstant);
514 
516 
518 
519 
521 
535  int ExtractCopy(int *A, int MyLDA) const;
536 
538 
548  int ExtractCopy(int **ArrayOfPointers) const;
549 
550  // ExtractView functions
551 
552 
554 
568  int ExtractView(int **A, int *MyLDA) const;
569 
571 
581  int ExtractView(int ***ArrayOfPointers) const;
582 
584 
586 
587 
589 
598  int MinValue (int * Result) const;
599 
601 
610  int MaxValue (int * Result) const;
611 
613 
615 
616 
618 
625 
626  // Local element access functions
627 
628  //
629 
631 
634  int*& operator [] (int i) { return Pointers_[i]; }
636 
639  // const int*& operator [] (int i) const;
640  int * const & operator [] (int i) const { return Pointers_[i]; }
641 
643 
646  Epetra_IntVector * & operator () (int i);
648 
651  const Epetra_IntVector * & operator () (int i) const;
652 
654 
656 
657 
659  int NumVectors() const {return(NumVectors_);}
660 
662  int MyLength() const {return(MyLength_);}
663 
665 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
666  int GlobalLength() const {
667  if(Map().GlobalIndicesInt())
668  return (int) GlobalLength_;
669  throw "Epetra_MultiVector::GlobalLength: GlobalIndices not int.";
670  }
671 #endif
672  long long GlobalLength64() const {return(GlobalLength_);}
673 
675  int Stride() const {return(Stride_);}
676 
678  bool ConstantStride() const {return(ConstantStride_);}
680 
684  int ReplaceMap(const Epetra_BlockMap& map);
685 
687 
688 
690  virtual void Print(std::ostream & os) const;
692 
694 
695 
697 
716  int ResetView(int ** ArrayOfPointers);
717 
719  int* Values() const {return Values_;}
720 
722  int** Pointers() const {return Pointers_;}
724 
725  // Expert-only function
726  // SuperLU defines Reduce to be a macro in util.h
727 #ifdef Reduce
728 #undef Reduce
729 #endif
730  int Reduce();
731 
732  protected:
733 
734  // Internal utilities
735  void Assign(const Epetra_IntMultiVector& rhs);
736  int CheckInput();
737 
738  int *Values_; // local MultiVector coefficients
739 
740  private:
741 
742 
743  // Internal utilities
744 
745  int AllocateForCopy(void);
746  int DoCopy(void);
747 
748  inline void UpdateOrdinalTemp() const
749  {if (OrdinalTemp_==0) OrdinalTemp_=new int[NumVectors_+1]; return;}
750 
751  inline void UpdateIntVectors() const {if (IntVectors_==0) { IntVectors_ = new Epetra_IntVector *[NumVectors_];
752  for (int i=0; i<NumVectors_; i++) IntVectors_[i] = 0;}
753  return;
754  }
755 
756  int AllocateForView(void);
757  int DoView(void);
758  template<typename int_type>
759  int ChangeGlobalValue(int_type GlobalBlockRow,
760  int BlockRowOffset,
761  int VectorIndex,
762  int OrdinalValue,
763  bool SumInto);
764  int ChangeMyValue(int MyBlockRow,
765  int BlockRowOffset,
766  int VectorIndex,
767  int OrdinalValue,
768  bool SumInto);
769 
770  int CheckSizes(const Epetra_SrcDistObject& A);
771 
772  int CopyAndPermute(const Epetra_SrcDistObject & Source,
773  int NumSameIDs,
774  int NumPermuteIDs,
775  int * PermuteToLIDs,
776  int * PermuteFromLIDs,
777  const Epetra_OffsetIndex * Indexor,
778  Epetra_CombineMode CombineMode = Zero);
779 
780  int PackAndPrepare(const Epetra_SrcDistObject & Source,
781  int NumExportIDs,
782  int * ExportLIDs,
783  int & LenExports,
784  char * & Exports,
785  int & SizeOfPacket,
786  int * Sizes,
787  bool & VarSizes,
788  Epetra_Distributor & Distor);
789 
790  int UnpackAndCombine(const Epetra_SrcDistObject & Source,
791  int NumImportIDs,
792  int * ImportLIDs,
793  int LenImports,
794  char * Imports,
795  int & SizeOfPacket,
796  Epetra_Distributor & Distor,
797  Epetra_CombineMode CombineMode,
798  const Epetra_OffsetIndex * Indexor );
799 
800  int **Pointers_; // Pointers to each vector;
801 
803  long long GlobalLength_;
807  int Stride_;
809  mutable int * OrdinalTemp_;
812 
813 };
814 
815 #endif /* EPETRA_INTMULTIVECTOR_H */
int ** Pointers() const
Get pointer to individual vector pointers.
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:127
Epetra_Distributor: The Epetra Gather/Scatter Setup Base Class.
virtual int CopyAndPermute(const Epetra_SrcDistObject &Source, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode=Zero)=0
Perform ID copies and permutations that are on processor.
Epetra_OffsetIndex: This class builds index for efficient mapping of data from one Epetra_CrsGraph ba...
Epetra_IntVector: A class for constructing and using dense integer vectors on a parallel computer...
int MyLength() const
Returns the local vector length on the calling processor of vectors in the multi-vector.
virtual void Print(std::ostream &os) const
Print method.
Epetra_Export: This class builds an export object for efficient exporting of off-processor elements...
Definition: Epetra_Export.h:70
bool ConstantStride() const
Returns true if this multi-vector has constant stride between vectors.
int NumVectors() const
Returns the number of vectors in the multi-vector.
int * Values() const
Get pointer to MultiVector values.
Epetra_IntMultiVector: A class for constructing and using dense multi-vectors, vectors and matrices i...
Epetra_CompObject & operator=(const Epetra_CompObject &src)
Epetra_Import: This class builds an import object for efficient importing of off-processor elements...
Definition: Epetra_Import.h:71
virtual int CheckSizes(const Epetra_SrcDistObject &Source)=0
Allows the source and target (this) objects to be compared for compatibility, return nonzero if not...
Epetra_IntVector ** IntVectors_
Epetra_Util: The Epetra Util Wrapper Class.
Definition: Epetra_Util.h:87
Epetra_BLAS: The Epetra BLAS Wrapper Class.
Definition: Epetra_BLAS.h:78
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:81
Epetra_CompObject: Functionality and data that is common to all computational classes.
virtual int UnpackAndCombine(const Epetra_SrcDistObject &Source, int NumImportIDs, int *ImportLIDs, int LenImports, char *Imports, int &SizeOfPacket, Epetra_Distributor &Distor, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor)=0
Perform any unpacking and combining after call to DoTransfer().
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
int GlobalLength() const
Returns the global vector length of vectors in the multi-vector.
int Stride() const
Returns the stride between vectors in the multi-vector (only meaningful if ConstantStride() is true)...
Epetra_CombineMode
Epetra_SrcDistObject: A class for supporting flexible source distributed objects for import/export op...
Epetra_DataAccess
long long GlobalLength64() const
Epetra_DistObject: A class for constructing and using dense multi-vectors, vectors and matrices in pa...
virtual int PackAndPrepare(const Epetra_SrcDistObject &Source, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor)=0
Perform any packing or preparation required for call to DoTransfer().
const Epetra_BlockMap & Map() const
Returns the address of the Epetra_BlockMap for this multi-vector.