Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_BasicRowMatrix.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_BASICROWMATRIX_H
45 #define EPETRA_BASICROWMATRIX_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 #include "Epetra_ConfigDefs.h"
56 #include "Epetra_RowMatrix.h"
57 #include "Epetra_Object.h"
58 #include "Epetra_CompObject.h"
59 #include "Epetra_Map.h"
60 #include "Epetra_Comm.h"
63 #include "Epetra_MultiVector.h"
64 
65 class Epetra_Vector;
66 class Epetra_Import;
67 class Epetra_Export;
68 
70 
118 class EPETRA_LIB_DLL_EXPORT Epetra_BasicRowMatrix: public Epetra_CompObject, public Epetra_Object, public virtual Epetra_RowMatrix {
119 
120  public:
121 
123 
124 
131  Epetra_BasicRowMatrix(const Epetra_Comm & Comm);
132 
134  virtual ~Epetra_BasicRowMatrix();
136 
138 
139 
151  void SetMaps(const Epetra_Map & RowMap, const Epetra_Map & ColMap);
152 
165  void SetMaps(const Epetra_Map & RowMap, const Epetra_Map & ColMap,
166  const Epetra_Map & DomainMap, const Epetra_Map & RangeMap);
167 
169 
170 
172 
173 
175 
184  virtual int ExtractMyRowCopy(int MyRow, int Length, int & NumEntries, double *Values, int * Indices) const = 0;
185 
187 
195  virtual int ExtractMyEntryView(int CurEntry, double * & Value, int & RowIndex, int & ColIndex) = 0;
196 
198 
206  virtual int ExtractMyEntryView(int CurEntry, double const * & Value, int & RowIndex, int & ColIndex) const = 0;
207 
209 
216  virtual int NumMyRowEntries(int MyRow, int & NumEntries) const = 0;
218 
220 
221 
223 
230  virtual int Multiply(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
231 
233 
242  virtual int Solve(bool Upper, bool Trans, bool UnitDiagonal,
243  const Epetra_MultiVector& X,
244  Epetra_MultiVector& Y) const
245  {
246  (void)Upper;
247  (void)Trans;
248  (void)UnitDiagonal;
249  (void)X;
250  (void)Y;
251  return(-1);
252  }
253 
255 
260  virtual int ExtractDiagonalCopy(Epetra_Vector & Diagonal) const;
261 
263 
272  virtual int InvRowSums(Epetra_Vector& x) const;
273 
275 
281  virtual int LeftScale(const Epetra_Vector& x);
282 
284 
293  virtual int InvColSums(Epetra_Vector& x) const;
294 
296 
302  virtual int RightScale(const Epetra_Vector& x);
304 
306 
307 
308 
310  virtual bool Filled() const {return(true);}
311 
313  bool LowerTriangular() const {if (!HaveNumericConstants_) ComputeNumericConstants(); return(LowerTriangular_);}
314 
316  virtual bool UpperTriangular() const {if (!HaveNumericConstants_) ComputeNumericConstants(); return(UpperTriangular_);}
317 
319 
321 
322 
331  virtual double NormInf() const{if (!HaveNumericConstants_) ComputeNumericConstants(); return(NormInf_);}
332 
341  virtual double NormOne() const{if (!HaveNumericConstants_) ComputeNumericConstants(); return(NormOne_);}
342 
349 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
350  virtual int NumGlobalNonzeros() const {
351  if(OperatorRangeMap().GlobalIndicesInt() && OperatorDomainMap().GlobalIndicesInt()) {
352  if (!HaveStructureConstants_)
353  ComputeStructureConstants();
354  return (int) NumGlobalNonzeros_;
355  }
356 
357  throw "Epetra_BasicRowMatrix::NumGlobalNonzeros: GlobalIndices not int.";
358  }
359 #endif
360  virtual long long NumGlobalNonzeros64() const{if (!HaveStructureConstants_) ComputeStructureConstants(); return(NumGlobalNonzeros_);}
361 
363 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
364  virtual int NumGlobalRows() const {
365  if(OperatorRangeMap().GlobalIndicesInt()) {
366  return (int) NumGlobalRows64();
367  }
368 
369  throw "Epetra_BasicRowMatrix::NumGlobalRows: GlobalIndices not int.";
370  }
371 #endif
372  virtual long long NumGlobalRows64() const {return(OperatorRangeMap().NumGlobalPoints64());}
373 
375 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
376  virtual int NumGlobalCols() const {
377  if(OperatorDomainMap().GlobalIndicesInt()) {
378  return (int) NumGlobalCols64();
379  }
380 
381  throw "Epetra_BasicRowMatrix::NumGlobalCols: GlobalIndices not int.";
382  }
383 #endif
384  virtual long long NumGlobalCols64() const {return(OperatorDomainMap().NumGlobalPoints64());}
385 
387 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
388  virtual int NumGlobalDiagonals() const {
389  if(OperatorDomainMap().GlobalIndicesInt()) {
390  return (int) NumGlobalDiagonals64();
391  }
392 
393  throw "Epetra_BasicRowMatrix::NumGlobalDiagonals: GlobalIndices not int.";
394  }
395 #endif
396  virtual long long NumGlobalDiagonals64() const{return(OperatorDomainMap().NumGlobalPoints64());}
397 
399  virtual int NumMyNonzeros() const{if (!HaveStructureConstants_) ComputeStructureConstants(); return(NumMyNonzeros_);}
400 
402  virtual int NumMyRows() const {return(OperatorRangeMap().NumMyPoints());}
403 
405  virtual int NumMyCols() const {return(RowMatrixColMap().NumMyPoints());}
406 
408  virtual int NumMyDiagonals() const {return(OperatorRangeMap().NumMyPoints());}
409 
411  virtual int MaxNumEntries() const{ if (!HaveStructureConstants_) ComputeStructureConstants(); return(MaxNumEntries_);}
412 
414  virtual const Epetra_Map & OperatorDomainMap() const {return(OperatorDomainMap_);}
415 
417  virtual const Epetra_Map & OperatorRangeMap() const {return(OperatorRangeMap_);}
418 
420  virtual const Epetra_BlockMap& Map() const {return(RowMatrixRowMap());}
421 
423  virtual const Epetra_Map & RowMatrixRowMap() const {return(RowMatrixRowMap_);}
424 
426  virtual const Epetra_Map & RowMatrixColMap() const {return(RowMatrixColMap_);}
427 
429  virtual const Epetra_Import * RowMatrixImporter() const {return(Importer_);}
430 
432  virtual const Epetra_Comm & Comm() const {return(*Comm_);}
434 
435 
437 
438 
440  virtual void Print(std::ostream & os) const;
442 
444 
445 
447 
455  virtual int SetUseTranspose(bool use_transpose) {UseTranspose_ = use_transpose; return(0);}
456 
458  virtual const char* Label() const {return(Epetra_Object::Label());}
459 
461 
467  virtual int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
469 
471 
479  virtual int ApplyInverse(const Epetra_MultiVector& X,
480  Epetra_MultiVector& Y) const
481  {
482  (void)X;
483  (void)Y;
484  return(-1);
485  }
486 
488  bool HasNormInf() const {return(true);}
489 
491  virtual bool UseTranspose() const {return(UseTranspose_);}
492 
494 
496 
497 
499 
507  virtual const Epetra_Import* Importer() const {return(Importer_);}
508 
510 
519  virtual const Epetra_Export* Exporter() const {return(Exporter_);}
520 
522 
523  protected:
524 
526 
527  /* Several constants are pre-computed to save excess computations. However, if the structure of the
529  problem changes, specifically if the nonzero count in any given row changes, then this function should be called
530  to update these constants.
531  */
532  virtual void ComputeStructureConstants() const;
534  /* Several numeric constants are pre-computed to save excess computations. However, if the values of the
535  problem change, then this function should be called to update these constants.
536  */
537  virtual void ComputeNumericConstants() const;
539 
540  void Setup();
541  void UpdateImportVector(int NumVectors) const;
542  void UpdateExportVector(int NumVectors) const;
543  void SetImportExport();
549 
550  mutable int NumMyNonzeros_;
551  mutable long long NumGlobalNonzeros_;
552  mutable int MaxNumEntries_;
553  mutable double NormInf_;
554  mutable double NormOne_;
557 
560  mutable bool LowerTriangular_;
561  mutable bool UpperTriangular_;
563  mutable bool HaveNumericConstants_;
564  mutable bool HaveMaps_;
565 
566 
571 
572 };
573 #endif /* EPETRA_BASICROWMATRIX_H */
Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in p...
virtual const Epetra_Map & RowMatrixRowMap() const
Returns the Row Map object needed for implementing Epetra_RowMatrix.
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:127
virtual const Epetra_BlockMap & Map() const
Implement the Epetra_SrcDistObjec::Map() function.
virtual long long NumGlobalDiagonals64() const =0
virtual int SetUseTranspose(bool use_transpose)
If set true, transpose of this operator will be applied.
virtual void Print(std::ostream &os) const
Print object to an output stream Print method.
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
virtual int RightScale(const Epetra_Vector &x)=0
Scales the Epetra_RowMatrix on the right with a Epetra_Vector x.
virtual double NormOne() const
Returns the one norm of the global matrix.
virtual long long NumGlobalDiagonals64() const
virtual const Epetra_Map & RowMatrixRowMap() const =0
Returns the Epetra_Map object associated with the rows of this matrix.
virtual int NumMyDiagonals() const
Returns the number of local nonzero diagonal entries.
virtual const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this matrix.
virtual long long NumGlobalNonzeros64() const
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
virtual bool Filled() const
If FillComplete() has been called, this query returns true, otherwise it returns false, presently always returns true.
virtual int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global matrix.
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_RowMatrix inverse applied to an Epetra_MultiVector X in Y...
Epetra_Export: This class builds an export object for efficient exporting of off-processor elements...
Definition: Epetra_Export.h:70
virtual const Epetra_Map & OperatorDomainMap() const =0
Returns the Epetra_Map object associated with the domain of this operator.
Epetra_Vector: A class for constructing and using dense vectors on a parallel computer.
virtual int ExtractDiagonalCopy(Epetra_Vector &Diagonal) const =0
Returns a copy of the main diagonal in a user-provided vector.
virtual int LeftScale(const Epetra_Vector &x)=0
Scales the Epetra_RowMatrix on the left with a Epetra_Vector x.
Epetra_Import: This class builds an import object for efficient importing of off-processor elements...
Definition: Epetra_Import.h:71
virtual int InvRowSums(Epetra_Vector &x) const =0
Computes the sum of absolute values of the rows of the Epetra_RowMatrix, results returned in x...
virtual int NumMyRowEntries(int MyRow, int &NumEntries) const =0
Returns the number of nonzero entries in MyRow.
virtual const Epetra_Import * RowMatrixImporter() const
Returns the Epetra_Import object that contains the import operations for distributed operations...
virtual const Epetra_Map & OperatorRangeMap() const =0
Returns the Epetra_Map object associated with the range of this operator.
virtual const char * Label() const
Epetra_Object Label access funtion.
virtual int MaxNumEntries() const
Returns the maximum number of nonzero entries across all rows on this processor.
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:81
virtual double NormInf() const
Returns the infinity norm of the global matrix.
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_RowMatrix applied to a Epetra_MultiVector X in Y.
Epetra_Object: The base Epetra class.
Definition: Epetra_Object.h:65
virtual int NumMyNonzeros() const
Returns the number of nonzero entries in the calling processor's portion of the matrix.
virtual bool UpperTriangular() const
If matrix is upper triangular, this query returns true, otherwise it returns false.
Epetra_MultiVector * ExportVector_
Epetra_CompObject: Functionality and data that is common to all computational classes.
virtual int NumGlobalRows() const
Returns the number of global matrix rows.
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
bool LowerTriangular() const
If matrix is lower triangular, this query returns true, otherwise it returns false.
Epetra_BasicRowMatrix: A class for simplifying the development of Epetra_RowMatrix adapters...
virtual const Epetra_Map & RowMatrixColMap() const
Returns the Column Map object needed for implementing Epetra_RowMatrix.
virtual long long NumGlobalCols64() const =0
virtual int NumGlobalCols() const
Returns the number of global matrix columns.
virtual int Solve(bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_BasicRowMatrix solve with a Epetra_MultiVector X in Y (not implemented...
virtual int Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y.
virtual const Epetra_Import * Importer() const
Returns the Epetra_Import object that contains the import operations for distributed operations...
virtual const char * Label() const
Returns a character string describing the operator.
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator (same as domain)...
virtual int NumMyRows() const
Returns the number of matrix rows owned by the calling processor.
virtual int NumMyCols() const
Returns the number of matrix columns owned by the calling processor.
virtual const Epetra_Export * Exporter() const
Returns the Epetra_Export object that contains the export operations for distributed operations...
virtual int ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const =0
Returns a copy of the specified local row in user-provided arrays.
virtual long long NumGlobalCols64() const
virtual long long NumGlobalRows64() const
virtual const Epetra_Map & RowMatrixColMap() const =0
Returns the Epetra_Map object associated with the columns of this matrix.
bool HasNormInf() const
Returns true because this class can compute an Inf-norm.
virtual int InvColSums(Epetra_Vector &x) const =0
Computes the sum of absolute values of the columns of the Epetra_RowMatrix, results returned in x...
Epetra_RowMatrix: A pure virtual class for using real-valued double-precision row matrices...
virtual int NumGlobalDiagonals() const
Returns the number of global nonzero diagonal entries.
virtual long long NumGlobalRows64() const =0
Epetra_MultiVector * ImportVector_