Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_LongLongSerialDenseMatrix.cpp
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 #include "Epetra_ConfigDefs.h"
45 #include "Epetra_Util.h"
46 
47 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
48 
49 //=============================================================================
51  : Epetra_Object("Epetra::LongLongSerialDenseMatrix"),
52  CV_(Copy),
53  A_Copied_(false),
54  M_(0),
55  N_(0),
56  LDA_(0),
57  A_(0)
58 {
59 }
60 
61 //=============================================================================
63  : Epetra_Object("Epetra::LongLongSerialDenseMatrix"),
64  CV_(Copy),
65  A_Copied_(false),
66  M_(0),
67  N_(0),
68  LDA_(0),
69  A_(0)
70 {
71  if(NumRows < 0)
72  throw ReportError("NumRows = " + toString(NumRows) + ". Should be >= 0", -1);
73  if(NumCols < 0)
74  throw ReportError("NumCols = " + toString(NumCols) + ". Should be >= 0", -1);
75 
76  int errorcode = Shape(NumRows, NumCols);
77  if(errorcode != 0)
78  throw ReportError("Shape returned non-zero (" + toString(errorcode) + ").", -2);
79 }
80 //=============================================================================
82  int NumRows, int NumCols)
83  : Epetra_Object("Epetra::LongLongSerialDenseMatrix"),
84  CV_(CV_in),
85  A_Copied_(false),
86  M_(NumRows),
87  N_(NumCols),
88  LDA_(lda),
89  A_(A_in)
90 {
91  if(A_in == 0)
92  throw ReportError("Null pointer passed as A_in parameter.", -3);
93  if(NumRows < 0)
94  throw ReportError("NumRows = " + toString(NumRows) + ". Should be >= 0", -1);
95  if(NumCols < 0)
96  throw ReportError("NumCols = " + toString(NumCols) + ". Should be >= 0", -1);
97  if(lda < 0)
98  throw ReportError("LDA = " + toString(lda) + ". Should be >= 0", -1);
99 
100  if(CV_in == Copy) {
101  LDA_ = M_;
102  const int newsize = LDA_ * N_;
103  if(newsize > 0) {
104  A_ = new long long[newsize];
105  CopyMat(A_in, lda, M_, N_, A_, LDA_);
106  A_Copied_ = true;
107  }
108  else {
109  A_ = 0;
110  }
111  }
112 }
113 //=============================================================================
115  : Epetra_Object(Source),
116  CV_(Source.CV_),
117  A_Copied_(false),
118  M_(Source.M_),
119  N_(Source.N_),
120  LDA_(Source.LDA_),
121  A_(Source.A_)
122 {
123  if(CV_ == Copy) {
124  LDA_ = M_;
125  const int newsize = LDA_ * N_;
126  if(newsize > 0) {
127  A_ = new long long[newsize];
128  CopyMat(Source.A_, Source.LDA_, M_, N_, A_, LDA_);
129  A_Copied_ = true;
130  }
131  else {
132  A_ = 0;
133  A_Copied_ = false;
134  }
135  }
136 }
137 //=============================================================================
138 int Epetra_LongLongSerialDenseMatrix::Reshape(int NumRows, int NumCols) {
139  if(NumRows < 0 || NumCols < 0)
140  return(-1);
141 
142  long long* A_tmp = 0;
143  const int newsize = NumRows * NumCols;
144 
145  if(newsize > 0) {
146  // Allocate space for new matrix
147  A_tmp = new long long[newsize];
148  for(int k = 0; k < newsize; k++)
149  A_tmp[k] = 0; // Zero out values
150  int M_tmp = EPETRA_MIN(M_, NumRows);
151  int N_tmp = EPETRA_MIN(N_, NumCols);
152  if(A_ != 0)
153  CopyMat(A_, LDA_, M_tmp, N_tmp, A_tmp, NumRows); // Copy principal submatrix of A to new A
154  }
155 
156  CleanupData(); // Get rid of anything that might be already allocated
157  M_ = NumRows;
158  N_ = NumCols;
159  LDA_ = M_;
160  A_ = A_tmp; // Set pointer to new A
161  A_Copied_ = (newsize>0);
162  return(0);
163 }
164 //=============================================================================
165 int Epetra_LongLongSerialDenseMatrix::Shape(int NumRows, int NumCols) {
166  if(NumRows < 0 || NumCols < 0)
167  return(-1);
168 
169  CleanupData(); // Get rid of anything that might be already allocated
170  M_ = NumRows;
171  N_ = NumCols;
172  LDA_ = M_;
173  const int newsize = LDA_ * N_;
174  if(newsize > 0) {
175  A_ = new long long[newsize];
176 #ifdef EPETRA_HAVE_OMP
177 #pragma omp parallel for
178 #endif
179  for(int k = 0; k < newsize; k++)
180  A_[k] = 0; // Zero out values
181  A_Copied_ = true;
182  }
183 
184  return(0);
185 }
186 //=============================================================================
188 {
189  CleanupData();
190 }
191 //=============================================================================
193 {
194  if(A_Copied_)
195  delete[] A_;
196  A_ = 0;
197  A_Copied_ = false;
198  M_ = 0;
199  N_ = 0;
200  LDA_ = 0;
201 }
202 //=============================================================================
204  if(this == &Source)
205  return(*this); // Special case of source same as target
206  if((CV_ == View) && (Source.CV_ == View) && (A_ == Source.A_))
207  return(*this); // Special case of both are views to same data.
208 
209  if(std::strcmp(Label(), Source.Label()) != 0)
210  throw ReportError("operator= type mismatch (lhs = " + std::string(Label()) +
211  ", rhs = " + std::string(Source.Label()) + ").", -5);
212 
213  if(Source.CV_ == View) {
214  if(CV_ == Copy) { // C->V only
215  CleanupData();
216  CV_ = View;
217  }
218  M_ = Source.M_; // C->V and V->V
219  N_ = Source.N_;
220  LDA_ = Source.LDA_;
221  A_ = Source.A_;
222  }
223  else {
224  if(CV_ == View) { // V->C
225  CV_ = Copy;
226  M_ = Source.M_;
227  N_ = Source.N_;
228  LDA_ = Source.M_;
229  const int newsize = LDA_ * N_;
230  if(newsize > 0) {
231  A_ = new long long[newsize];
232  A_Copied_ = true;
233  }
234  else {
235  A_ = 0;
236  A_Copied_ = false;
237  }
238  }
239  else { // C->C
240  if((Source.M_ <= LDA_) && (Source.N_ == N_)) { // we don't need to reallocate
241  M_ = Source.M_;
242  N_ = Source.N_;
243  }
244  else { // we need to allocate more space (or less space)
245  CleanupData();
246  M_ = Source.M_;
247  N_ = Source.N_;
248  LDA_ = Source.M_;
249  const int newsize = LDA_ * N_;
250  if(newsize > 0) {
251  A_ = new long long[newsize];
252  A_Copied_ = true;
253  }
254  }
255  }
256  CopyMat(Source.A_, Source.LDA_, M_, N_, A_, LDA_); // V->C and C->C
257  }
258 
259  return(*this);
260 }
261 
262 //=============================================================================
264 {
265  if (M_ != rhs.M_ || N_ != rhs.N_) return(false);
266 
267  const long long* A_tmp = A_;
268  const long long* rhsA = rhs.A_;
269 
270  for(int j=0; j<N_; ++j) {
271  int offset = j*LDA_;
272  int rhsOffset = j*rhs.LDA_;
273  for(int i=0; i<M_; ++i) {
274  if (A_tmp[offset+i] != rhsA[rhsOffset+i]) {
275  return(false);
276  }
277  }
278  }
279 
280  return(true);
281 }
282 
283 //=============================================================================
285  if(std::strcmp(Label(), Source.Label()) != 0)
286  return(-1);
287 
288  if(CV_ == Copy) {
289  CleanupData();
290  CV_ = View;
291  }
292  M_ = Source.M_;
293  N_ = Source.N_;
294  LDA_ = Source.LDA_;
295  A_ = Source.A_;
296 
297  return(0);
298 }
299 
300 //=============================================================================
301 void Epetra_LongLongSerialDenseMatrix::CopyMat(long long* Source, int Source_LDA, int NumRows, int NumCols,
302  long long* Target, int Target_LDA)
303 {
304  int i, j;
305  long long* targetPtr = Target;
306  long long* sourcePtr = 0;
307  for(j = 0; j < NumCols; j++) { // for each column
308  targetPtr = Target + j * Target_LDA; // set pointers to correct stride
309  sourcePtr = Source + j * Source_LDA;
310  for(i = 0; i < NumRows; i++) // for each row
311  *targetPtr++ = *sourcePtr++; // copy element (i,j) and increment pointer to (i,j+1)
312  }
313  return;
314 }
315 
316 //=============================================================================
318  long long anorm = 0;
319  long long* ptr = 0;
320  for(int j = 0; j < N_; j++) {
321  long long sum = 0;
322  ptr = A_ + j*LDA_;
323  for(int i = 0; i < M_; i++)
324  {
325  const long long val = *ptr++;
326  sum += (val > 0 ? val : -val); // No std::abs(long long) on VS2005.
327  }
328  anorm = EPETRA_MAX(anorm, sum);
329  }
330  return(anorm);
331 }
332 
333 //=============================================================================
335  long long anorm = 0;
336  long long* ptr = 0;
337  // Loop across columns in inner loop. Most expensive memory access, but
338  // requires no extra storage.
339  for(int i = 0; i < M_; i++) {
340  long long sum = 0;
341  ptr = A_ + i;
342  for(int j = 0; j < N_; j++) {
343  const long long val = *ptr;
344  sum += (val > 0 ? val : -val); // No std::abs(long long) on VS2005.
345  ptr += LDA_;
346  }
347  anorm = EPETRA_MAX(anorm, sum);
348  }
349  return(anorm);
350 }
351 
352 //=========================================================================
353 void Epetra_LongLongSerialDenseMatrix::Print(std::ostream& os) const {
354  if(CV_ == Copy)
355  os << "Data access mode: Copy" << std::endl;
356  else
357  os << "Data access mode: View" << std::endl;
358  if(A_Copied_)
359  os << "A_Copied: yes" << std::endl;
360  else
361  os << "A_Copied: no" << std::endl;
362  os << "Rows(M): " << M_ << std::endl;
363  os << "Columns(N): " << N_ << std::endl;
364  os << "LDA: " << LDA_ << std::endl;
365  if(M_ == 0 || N_ == 0)
366  os << "(matrix is empty, no values to display)" << std::endl;
367  else
368  for(int i = 0; i < M_; i++) {
369  for(int j = 0; j < N_; j++){
370  os << (*this)(i,j) << " ";
371  }
372  os << std::endl;
373  }
374 }
375 
376 //=========================================================================
378 
379  Epetra_Util util;
380 
381  for(int j = 0; j < N_; j++) {
382  long long* arrayPtr = A_ + (j * LDA_);
383  for(int i = 0; i < M_; i++) {
384  *arrayPtr++ = util.RandomInt();
385  }
386  }
387 
388  return(0);
389 }
390 
391 #endif // EPETRA_NO_64BIT_GLOBAL_INDICES
Epetra_LongLongSerialDenseMatrix()
Default constructor; defines a zero size object.
int MakeViewOf(const Epetra_LongLongSerialDenseMatrix &Source)
Reset an existing LongLongSerialDenseMatrix to point to another Matrix.
bool operator==(const Epetra_LongLongSerialDenseMatrix &rhs) const
Comparison operator.
virtual long long InfNorm()
Computes the Infinity-Norm of the this matrix.
#define EPETRA_MIN(x, y)
void CopyMat(long long *Source, int Source_LDA, int NumRows, int NumCols, long long *Target, int Target_LDA)
unsigned int RandomInt()
Returns a random integer on the interval (0, 2^31-1)
Definition: Epetra_Util.cpp:71
virtual ~Epetra_LongLongSerialDenseMatrix()
Epetra_LongLongSerialDenseMatrix destructor.
virtual long long OneNorm()
Computes the 1-Norm of the this matrix.
Epetra_Util: The Epetra Util Wrapper Class.
Definition: Epetra_Util.h:87
virtual const char * Label() const
Epetra_Object Label access funtion.
std::string toString(const int &x) const
Epetra_Object: The base Epetra class.
Definition: Epetra_Object.h:65
int Random()
Set matrix values to random numbers.
int Shape(int NumRows, int NumCols)
Set dimensions of a Epetra_LongLongSerialDenseMatrix object; init values to zero. ...
int Reshape(int NumRows, int NumCols)
Reshape a Epetra_LongLongSerialDenseMatrix object.
Epetra_LongLongSerialDenseMatrix & operator=(const Epetra_LongLongSerialDenseMatrix &Source)
Copy from one matrix to another.
virtual int ReportError(const std::string Message, int ErrorCode) const
Error reporting method.
virtual void Print(std::ostream &os) const
Print service methods; defines behavior of ostream &lt;&lt; operator.
Epetra_DataAccess
#define EPETRA_MAX(x, y)
Epetra_LongLongSerialDenseMatrix: A class for constructing and using general dense integer matrices...