Epetra Package Browser (Single Doxygen Collection)  Development
Epetra_LongLongSerialDenseMatrix.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_LONGLONGSERIALDENSEMATRIX_H
45 #define EPETRA_LONGLONGSERIALDENSEMATRIX_H
46 
47 #include "Epetra_ConfigDefs.h"
48 #include "Epetra_Object.h"
49 
50 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
51 
53 
119 //=========================================================================
120 class EPETRA_LIB_DLL_EXPORT Epetra_LongLongSerialDenseMatrix : public Epetra_Object {
121 
122  public:
123 
125 
126 
133 
135 
146  Epetra_LongLongSerialDenseMatrix(int NumRows, int NumCols);
147 
149 
164  Epetra_LongLongSerialDenseMatrix(Epetra_DataAccess CV, long long* A, int LDA, int NumRows, int NumCols);
165 
167 
171 
175 
177 
178 
191  int Shape(int NumRows, int NumCols);
192 
194 
207  int Reshape(int NumRows, int NumCols);
209 
211 
212 
214 
217  virtual long long OneNorm();
218 
220  virtual long long InfNorm();
221 
223 
230 
232 
235  bool operator==(const Epetra_LongLongSerialDenseMatrix& rhs) const;
236 
238 
241  { return !(*this == rhs); }
242 
244 
253  long long& operator () (int RowIndex, int ColIndex);
254 
256 
265  const long long& operator () (int RowIndex, int ColIndex) const;
266 
268 
278  long long* operator [] (int ColIndex);
279 
281 
291  const long long* operator [] (int ColIndex) const;
292 
294 
300  int Random();
301 
303  int M() const {return(M_);};
304 
306  int N() const {return(N_);};
307 
309  const long long* A() const {return(A_);};
310 
312  long long* A() {return(A_);};
313 
315  int LDA() const {return(LDA_);};
316 
318  Epetra_DataAccess CV() const {return(CV_);};
320 
322 
323  virtual void Print(std::ostream& os) const;
326 
328 
329 
331 
348  int MakeViewOf(const Epetra_LongLongSerialDenseMatrix& Source);
350 
351  protected:
352 
353  void CopyMat(long long* Source, int Source_LDA, int NumRows, int NumCols, long long* Target, int Target_LDA);
354  void CleanupData();
355 
357  bool A_Copied_;
358  int M_;
359  int N_;
360  int LDA_;
361  long long* A_;
362 
363 };
364 
365 // inlined definitions of op() and op[]
366 //=========================================================================
367 inline long long& Epetra_LongLongSerialDenseMatrix::operator () (int RowIndex, int ColIndex) {
368 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
369  if(RowIndex >= M_ || RowIndex < 0)
370  throw ReportError("Row index = " + toString(RowIndex) +
371  " Out of Range 0 - " + toString(M_-1),-1);
372  if(ColIndex >= N_ || ColIndex < 0)
373  throw ReportError("Column index = " + toString(ColIndex) +
374  " Out of Range 0 - " + toString(N_-1),-2);
375 #endif
376  return(A_[ColIndex*LDA_ + RowIndex]);
377 }
378 //=========================================================================
379 inline const long long& Epetra_LongLongSerialDenseMatrix::operator () (int RowIndex, int ColIndex) const {
380 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
381  if(RowIndex >= M_ || RowIndex < 0)
382  throw ReportError("Row index = " + toString(RowIndex) +
383  " Out of Range 0 - " + toString(M_-1),-1);
384  if(ColIndex >= N_ || ColIndex < 0)
385  throw ReportError("Column index = " + toString(ColIndex) +
386  " Out of Range 0 - " + toString(N_-1),-2);
387 #endif
388  return(A_[ColIndex * LDA_ + RowIndex]);
389 }
390 //=========================================================================
391 inline long long* Epetra_LongLongSerialDenseMatrix::operator [] (int ColIndex) {
392 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
393  if(ColIndex >= N_ || ColIndex < 0)
394  throw ReportError("Column index = " + toString(ColIndex) +
395  " Out of Range 0 - " + toString(N_-1),-2);
396 #endif
397  return(A_+ ColIndex * LDA_);
398 }
399 //=========================================================================
400 inline const long long* Epetra_LongLongSerialDenseMatrix::operator [] (int ColIndex) const {
401 #ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
402  if(ColIndex >= N_ || ColIndex < 0)
403  throw ReportError("Column index = " + toString(ColIndex) +
404  " Out of Range 0 - " + toString(N_-1),-2);
405 #endif
406  return(A_ + ColIndex * LDA_);
407 }
408 //=========================================================================
409 
410 #endif // EPETRA_NO_64BIT_GLOBAL_INDICES
411 
412 #endif /* EPETRA_LONGLONGSERIALDENSEMATRIX_H */
Epetra_DataAccess
Epetra_LongLongSerialDenseMatrix: A class for constructing and using general dense integer matrices.
int N() const
Returns column dimension of system.
int M() const
Returns row dimension of system.
int LDA() const
Returns the leading dimension of the this matrix.
bool operator!=(const Epetra_LongLongSerialDenseMatrix &rhs) const
Inequality operator.
long long * A()
Returns pointer to the this matrix.
Epetra_DataAccess CV() const
Returns the data access mode of the this matrix.
const long long * A() const
Returns const pointer to the this matrix.
long long * operator[](int ColIndex)
Column access function.
long long & operator()(int RowIndex, int ColIndex)
Element access function.
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.
virtual int ReportError(const std::string Message, int ErrorCode) const
Error reporting method.
std::string toString(const int &x) const
Epetra_Object & operator=(const Epetra_Object &src)