EpetraExt Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
EpetraExt_BlockDiagMatrix.h
Go to the documentation of this file.
1 //@HEADER
2 // ***********************************************************************
3 //
4 // EpetraExt: Epetra Extended - Linear Algebra Services Package
5 // Copyright (2011) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 
42 #ifndef EPETRAEXT_BLOCKDIAGMATRIX_H
43 #define EPETRAEXT_BLOCKDIAGMATRIX_H
44 
45 #include "Epetra_BLAS.h"
46 #include "Epetra_LAPACK.h"
47 #include "Epetra_DistObject.h"
48 #include "Epetra_BlockMap.h"
49 #include "Epetra_Map.h"
50 #include "Epetra_Operator.h"
52 
53 class Epetra_Comm;
54 
56 
65 //=========================================================================
67 
68  public:
70  EpetraExt_BlockDiagMatrix(const Epetra_BlockMap& Map,bool zero_out=true);
71 
72 
75 
78 
80 
87 
89 
92  double* operator [] (int index) {return &Values_[DataMap_->FirstPointInElement(index)];}
94 
97  const double* operator [] (int index) const {return &Values_[DataMap_->FirstPointInElement(index)];}
99 
100 
102 
103 
105  virtual int SetUseTranspose(bool /* useTranspose */){return -1;}
106 
108  virtual int SetParameters(Teuchos::ParameterList & List);
109 
111  virtual int Compute();
112 
114 
115 
117 
118 
120  virtual const char * Label() const{return "EpetraExt::BlockDiagMatrix";}//HAQ
121 
123  virtual bool UseTranspose() const {return false;}
124 
126  virtual bool HasNormInf() const {return false;}
127 
129  virtual const Epetra_Comm & Comm() const {return Map().Comm();}
130 
132  virtual const Epetra_Map & OperatorDomainMap() const {return *dynamic_cast<const Epetra_Map*>(&Map());}
133 
135  virtual const Epetra_Map & OperatorRangeMap() const {return *dynamic_cast<const Epetra_Map*>(&Map());}
136 
138  virtual const Epetra_BlockMap & BlockMap() const {return Map();}
139 
141  double* Values() const {return(Values_);}
142 
144  int BlockSize(int LID) const {return Map().ElementSize(LID);}
145 
147  int DataSize(int LID) const {return DataMap_->ElementSize(LID);}
148 
150  bool ConstantBlockSize() const {return Map().ConstantElementSize();}
151 
153  int NumMyBlocks() const {return(Map().NumMyElements());}
154 
156 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
157  int NumGlobalBlocks() const {return(Map().NumGlobalElements());}
158 #endif
159  long long NumGlobalBlocks64() const {return(Map().NumGlobalElements64());}
160 
162  int NumMyUnknowns() const {return(Map().NumMyPoints());}
163 
165 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
166  int NumGlobalUnknowns() const {return(Map().NumGlobalPoints());}
167 #endif
168  long long NumGlobalUnknowns64() const {return(Map().NumGlobalPoints64());}
169 
171  int NumData() const {return DataMap_->NumMyPoints();}
172 
174  int GetApplyMode() {return ApplyMode_;}
175 
177  virtual void Print(std::ostream & os) const;
178 
180 
181 
183 
184 
186 
193  virtual int Apply(const Epetra_MultiVector& /* X */, Epetra_MultiVector& /* Y */) const {return -1;}
194 
196 
207  virtual int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
208 
210  virtual double NormInf() const{return -1;}
211 
213  void PutScalar(double value);
214 
216  virtual const Epetra_BlockMap & DataMap() const {return *DataMap_;}
217 
219 
220 
221 private:
222  void Allocate();
223 
224  int DoCopy(const EpetraExt_BlockDiagMatrix& Source);
225 
226  // Routines to implement Epetra_DistObject virtual methods
227  // Allows the source and target (\e this) objects to be compared for compatibility, return nonzero if not.
228  int CheckSizes(const Epetra_SrcDistObject& Source);
229  // Perform ID copies and permutations that are on processor.
230  int CopyAndPermute(const Epetra_SrcDistObject& Source,
231  int NumSameIDs,
232  int NumPermuteIDs,
233  int * PermuteToLIDs,
234  int * PermuteFromLIDs,
235  const Epetra_OffsetIndex * Indexor,
236  Epetra_CombineMode CombineMode = Zero);
237 
238  // Perform any packing or preparation required for call to DoTransfer().
239  int PackAndPrepare(const Epetra_SrcDistObject& Source,
240  int NumExportIDs,
241  int* ExportLIDs,
242  int& LenExports,
243  char*& Exports,
244  int& SizeOfPacket,
245  int* Sizes,
246  bool & VarSizes,
247  Epetra_Distributor& Distor);
248 
249  // Perform any unpacking and combining after call to DoTransfer().
250  int UnpackAndCombine(const Epetra_SrcDistObject& Source,
251  int NumImportIDs,
252  int* ImportLIDs,
253  int LenImports,
254  char* Imports,
255  int& SizeOfPacket,
256  Epetra_Distributor& Distor,
257  Epetra_CombineMode CombineMode,
258  const Epetra_OffsetIndex * Indexor);
259 
262 
265 
268 
271 
273  double *Values_;
274 
276  int *Pivots_;
277 
278 }; /* EPETRAEXT_BLOCKDIAGMATRIX_H */
279 
280 #endif
virtual const Epetra_BlockMap & DataMap() const
Returns the Epetra_BlockMap object with the distribution of underlying values.
bool ConstantBlockSize() const
Returns true if the element size is constant.
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y...
int NumGlobalUnknowns() const
Returns the number of global unknowns.
int ElementSize() const
int * Pivots_
Pivots for factorization.
int NumMyBlocks() const
Returns the number of local blocks.
int PackAndPrepare(const Epetra_SrcDistObject &Source, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor)
Epetra_BlockMap * DataMap_
Map for the data.
double * Values() const
Returns a pointer to the array containing the blocks.
bool ConstantElementSize() const
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
virtual const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
int NumMyUnknowns() const
Returns the number of local unknowns.
EpetraExt_BlockDiagMatrix & operator=(const EpetraExt_BlockDiagMatrix &Source)
= Operator.
int NumGlobalBlocks() const
Returns the number of global blocks.
virtual double NormInf() const
NormInf - Not Implemented.
virtual ~EpetraExt_BlockDiagMatrix()
Destructor.
double * operator[](int index)
Block access function.
EpetraExt_BlockDiagMatrix: A class for storing distributed block matrices.
int BlockSize(int LID) const
Returns the size of the given block.
virtual void Print(std::ostream &os) const
Print method.
int NumData() const
Returns the size of the total Data block.
int FirstPointInElement(int LID) const
EpetraExt_BlockDiagMatrix(const Epetra_BlockMap &Map, bool zero_out=true)
Constructor - This map is the map of the vector this can be applied to.
virtual int SetUseTranspose(bool)
SetUseTranspose - not implemented.
virtual const Epetra_BlockMap & Map() const =0
int CopyAndPermute(const Epetra_SrcDistObject &Source, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode=Zero)
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)
virtual int Compute()
Computes the inverse / factorization if such is set on the list.
virtual int Apply(const Epetra_MultiVector &, Epetra_MultiVector &) const
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
virtual int SetParameters(Teuchos::ParameterList &List)
SetParameters.
const Epetra_Comm & Comm() const
double * Values_
Actual Data values.
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
bool HasComputed_
Has Computed? Needed for Inverse/Factorization modes.
int DoCopy(const EpetraExt_BlockDiagMatrix &Source)
virtual const char * Label() const
Returns a character std::string describing the operator.
Epetra_CombineMode
int ApplyMode_
Which Apply Mode to use.
int DataSize(int LID) const
Returns the size of the data in the given block.
int NumMyPoints() const
int GetApplyMode()
Gets apply mode info.
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
int CheckSizes(const Epetra_SrcDistObject &Source)
virtual const Epetra_BlockMap & BlockMap() const
Returns the Epetra_BlockMap object associated with the range of this operator.
void PutScalar(double value)
PutScalar function.