Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_SerialDenseMatrix.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_SERIALDENSEMATRIX_H
45 #define EPETRA_SERIALDENSEMATRIX_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_Object.h"
57 #include "Epetra_CompObject.h"
58 #include "Epetra_BLAS.h"
61 class Epetra_VbrMatrix;
62 
64 
114 //=========================================================================
115 class EPETRA_LIB_DLL_EXPORT Epetra_SerialDenseMatrix : public Epetra_CompObject, public Epetra_Object, public Epetra_SerialDenseOperator, public Epetra_BLAS {
116 
117  public:
118 
120 
121 
127  Epetra_SerialDenseMatrix(bool set_object_label=true);
128 
130 
141  Epetra_SerialDenseMatrix(int NumRows, int NumCols, bool set_object_label=true);
142 
144 
159  Epetra_SerialDenseMatrix(Epetra_DataAccess CV, double* A_in, int LDA_in, int NumRows, int NumCols,
160  bool set_object_label=true);
161 
163 
165 
167  virtual ~Epetra_SerialDenseMatrix ();
169 
171 
172 
185  int Shape(int NumRows, int NumCols);
186 
188 
201  int Reshape(int NumRows, int NumCols);
203 
205 
206 
208 
227  int Multiply(char TransA, char TransB, double ScalarAB,
230  double ScalarThis);
231 
233  /* This method is intended to imitate the semantics of the matrix-vector
234  multiplication provided by Epetra's sparse matrices. The 'vector' arguments
235  are actually matrices; this method will return an error if the
236  dimensions of 'x' are not compatible. 'y' will be reshaped if necessary.
237  */
238  int Multiply(bool transA,
239  const Epetra_SerialDenseMatrix& x,
241 
243 
264  int Multiply(char SideA, double ScalarAB,
266  const Epetra_SerialDenseMatrix& B,
267  double ScalarThis);
268 
270 
278  int Scale(double ScalarA);
279 
281 
284  virtual double NormOne() const;
285 
287  virtual double NormInf() const;
288 
290 
292 
293 
295 
302 
304 
307  bool operator==(const Epetra_SerialDenseMatrix& rhs) const;
308 
310 
312  bool operator!=(const Epetra_SerialDenseMatrix& rhs) const
313  { return !(*this == rhs); }
314 
316 
323  Epetra_SerialDenseMatrix & operator += (const Epetra_SerialDenseMatrix& Source);
324 
326 
335  double& operator () (int RowIndex, int ColIndex);
336 
338 
347  const double& operator () (int RowIndex, int ColIndex) const;
348 
350 
360  double* operator [] (int ColIndex);
361 
363 
373  const double* operator [] (int ColIndex) const;
374 
376 
382  int Random();
383 
385  int M() const {return(M_);};
386 
388  int N() const {return(N_);};
389 
391  double* A() const {return(A_);};
392 
394  double* A() {return(A_);};
395 
397  int LDA() const {return(LDA_);};
398 
400  Epetra_DataAccess CV() const {return(CV_);};
402 
404 
405  virtual void Print(std::ostream& os) const;
408 
410 
411 
413 
416  virtual double OneNorm() const {return(NormOne());};
417 
419  virtual double InfNorm() const {return(NormInf());};
421 
423 
424 
426 
435  virtual int SetUseTranspose(bool UseTranspose_in) { UseTranspose_ = UseTranspose_in; return (0); }
436 
438 
446  virtual int Apply(const Epetra_SerialDenseMatrix& X, Epetra_SerialDenseMatrix& Y);
447 
449 
459  {
460  (void)X;//prevents unused variable compiler warning
461  (void)Y;
462  return (-1);
463  }
464 
466  virtual const char * Label() const { return Epetra_Object::Label(); }
467 
469  virtual bool UseTranspose() const { return UseTranspose_; }
470 
472  virtual bool HasNormInf() const { return true; }
473 
475  virtual int RowDim() const { return M(); }
476 
478  virtual int ColDim() const { return N(); }
480 
481  protected:
482 
483  void CopyMat(const double* Source, int Source_LDA, int NumRows, int NumCols,
484  double* Target, int Target_LDA, bool add=false);
485  void CleanupData();
486 
487  int M_;
488  int N_;
489  bool A_Copied_;
491 
492  //For performance reasons, it's better if Epetra_VbrMatrix can access the
493  //LDA_ and A_ members of this class directly without going through an
494  //accessor method. Rather than making them public members, we'll make
495  //Epetra_VbrMatrix a friend class.
496 
497  friend class Epetra_VbrMatrix;
498 
499  int LDA_;
500  double* A_;
501 
503 };
504 
505 // inlined definitions of op() and op[]
506 //=========================================================================
507 inline double& Epetra_SerialDenseMatrix::operator () (int RowIndex, int ColIndex) {
508 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
509  if (RowIndex >= M_ || RowIndex < 0)
510  throw ReportError("Row index = " +toString(RowIndex) +
511  " Out of Range 0 - " + toString(M_-1),-1);
512  if (ColIndex >= N_ || ColIndex < 0)
513  throw ReportError("Column index = " +toString(ColIndex) +
514  " Out of Range 0 - " + toString(N_-1),-2);
515 #endif
516  return(A_[ColIndex*LDA_ + RowIndex]);
517 }
518 //=========================================================================
519 inline const double& Epetra_SerialDenseMatrix::operator () (int RowIndex, int ColIndex) const {
520 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
521  if (RowIndex >= M_ || RowIndex < 0)
522  throw ReportError("Row index = " +toString(RowIndex) +
523  " Out of Range 0 - " + toString(M_-1),-1);
524  if (ColIndex >= N_ || ColIndex < 0)
525  throw ReportError("Column index = " +toString(ColIndex) +
526  " Out of Range 0 - " + toString(N_-1),-2);
527 #endif
528  return(A_[ColIndex*LDA_ + RowIndex]);
529 }
530 //=========================================================================
531 inline double* Epetra_SerialDenseMatrix::operator [] (int ColIndex) {
532 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
533  if (ColIndex >= N_ || ColIndex < 0)
534  throw ReportError("Column index = " +toString(ColIndex) +
535  " Out of Range 0 - " + toString(N_-1),-2);
536 #endif
537  return(A_ + ColIndex*LDA_);
538 }
539 //=========================================================================
540 inline const double* Epetra_SerialDenseMatrix::operator [] (int ColIndex) const {
541 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
542  if (ColIndex >= N_ || ColIndex < 0)
543  throw ReportError("Column index = " +toString(ColIndex) +
544  " Out of Range 0 - " + toString(N_-1),-2);
545 #endif
546  return(A_+ ColIndex*LDA_);
547 }
548 //=========================================================================
549 
550 #endif /* EPETRA_SERIALDENSEMATRIX_H */
virtual void Print(std::ostream &os) const
Print object to an output stream Print method.
double * operator[](int ColIndex)
Column access function.
virtual int ColDim() const
Returns the column dimension of operator.
virtual int RowDim() const
Returns the row dimension of operator.
Epetra_DataAccess CV() const
Returns the data access mode of the this matrix.
Epetra_VbrMatrix: A class for the construction and use of real-valued double-precision variable block...
virtual int ApplyInverse(const Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &Y)
Returns the result of a Epetra_SerialDenseOperator inverse applied to an Epetra_SerialDenseMatrix X i...
Epetra_CompObject & operator=(const Epetra_CompObject &src)
virtual int Apply(const Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &Y)=0
Returns the result of a Epetra_SerialDenseOperator applied to a Epetra_SerialDenseMatrix X in Y...
Epetra_SerialDenseMatrix: A class for constructing and using real double precision general dense matr...
Epetra_SerialSymDenseMatrix: A class for constructing and using symmetric positive definite dense mat...
double * A() const
Returns pointer to the this matrix.
Epetra_BLAS: The Epetra BLAS Wrapper Class.
Definition: Epetra_BLAS.h:78
virtual const char * Label() const
Epetra_Object Label access funtion.
std::string toString(const int &x) const
double & operator()(int RowIndex, int ColIndex)
Element access function.
virtual double OneNorm() const
Computes the 1-Norm of the this matrix (identical to NormOne() method).
Epetra_Object: The base Epetra class.
Definition: Epetra_Object.h:65
virtual int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
Epetra_CompObject: Functionality and data that is common to all computational classes.
virtual double InfNorm() const
Computes the Infinity-Norm of the this matrix (identical to NormInf() method).
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Epetra_SerialDenseOperator: A pure virtual class for using real-valued double-precision operators...
virtual double NormInf() const =0
Returns the infinity norm of the global matrix.
int LDA() const
Returns the leading dimension of the this matrix.
bool operator!=(const Epetra_SerialDenseMatrix &rhs) const
Inequality operator.
virtual const char * Label() const
Returns a character string describing the operator.
double * A()
Returns pointer to the this matrix.
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
virtual int ReportError(const std::string Message, int ErrorCode) const
Error reporting method.
int N() const
Returns column dimension of system.
Epetra_DataAccess
bool operator==(BigUInt< n > const &a, BigUInt< n > const &b)
int M() const
Returns row dimension of system.