IFPACK Development
Loading...
Searching...
No Matches
Ifpack_SerialTriDiMatrix.h
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 IFPACK_SERIALTRIDIMATRIX_H
45#define IFPACK_SERIALTRIDIMATRIX_H
46
47#if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
48#ifdef __GNUC__
49#warning "The Ifpack package is deprecated"
50#endif
51#endif
52
53#include "Epetra_ConfigDefs.h"
54#include "Epetra_Object.h"
55#include "Epetra_CompObject.h"
56#include "Epetra_BLAS.h"
57
59
61
111//=========================================================================
113
114 public:
115
117
118
124 Ifpack_SerialTriDiMatrix(bool set_object_label=true);
125
127
136 Ifpack_SerialTriDiMatrix(int NumRowCol, bool set_object_label=true);
137
139
150 Ifpack_SerialTriDiMatrix(Epetra_DataAccess CV, double* A_in, int NumRowCol,
151 bool set_object_label=true);
152
154
156
158 virtual ~Ifpack_SerialTriDiMatrix ();
160
162
163
174 int Shape(int NumRowCol);
175
176 int Reshape(int, int);
177
179
180
182
201 int Multiply(char TransA, char TransB, double ScalarAB,
204 double ScalarThis);
205
207 /* This method is intended to imitate the semantics of the matrix-vector
208 multiplication provided by Epetra's sparse matrices. The 'vector' arguments
209 are actually matrices; this method will return an error if the
210 dimensions of 'x' are not compatible. 'y' will be reshaped if necessary.
211 */
212 /* int Multiply(bool transA, */
213 /* const Ifpack_SerialTriDiMatrix& x, */
214 /* Ifpack_SerialTriDiMatrix& y); */
215
217
244 int Scale(double ScalarA);
245
247
250 virtual double NormOne() const;
251
253 virtual double NormInf() const;
254
256
258
259
261
268
270
273 bool operator==(const Ifpack_SerialTriDiMatrix& rhs) const;
274
276
279 { return !(*this == rhs); }
280
282
290
292
300 double& operator () (int RowIndex, int ColIndex);
301
303
312 const double& operator () (int RowIndex, int ColIndex) const;
313
315
325 // double* operator [] (int ColIndex);
326
328
338 // const double* operator [] (int ColIndex) const;
339
341
347 int Random();
348
350 int N() const {return(N_);};
351
352 int LDA() const {return(LDA_);};
353
355 double* A() const {return(A_);};
356
358 // double* A() {return(A_);};
359
360 double* DL() { return DL_;};
361 double* DL() const { return DL_;};
362 double* D() { return D_;};
363 double* D() const { return D_;};
364 double* DU() { return DU_;};
365 double* DU() const { return DU_;};
366 double* DU2() { return DU2_;};
367 double* DU2() const { return DU2_;};
368
370 Epetra_DataAccess CV() const {return(CV_);};
372
374
375
376 virtual void Print(std::ostream& os) const;
378
380
381
383
386 virtual double OneNorm() const {return(NormOne());};
387
389 virtual double InfNorm() const {return(NormInf());};
391
393
394
396
405 virtual int SetUseTranspose(bool UseTranspose_in) { UseTranspose_ = UseTranspose_in; return (0); }
406
408
418 {
419 (void)X;//prevents unused variable compiler warning
420 (void)Y;
421 return (-1);
422 }
423
425 virtual const char * Label() const { return Epetra_Object::Label(); }
426
428 virtual bool UseTranspose() const { return UseTranspose_; }
429
431 virtual bool HasNormInf() const { return true; }
432
434 virtual int RowColDim() const { return N(); }
436
437 protected:
438
439 void CopyMat(const double* Source, int NumRowCol,
440 double* Target, int NRC2, bool add=false);
441 void CleanupData();
442
443 int N_;
444 int LDA_;
445 bool A_Copied_;
447
448 //For performance reasons, it's better if Epetra_VbrMatrix can access the
449 //A_ members of this class directly without going through an
450 //accessor method. Rather than making them public members, we'll make
451 //Epetra_VbrMatrix a friend class.
452
453 friend class Epetra_VbrMatrix;
454
455 double* A_;
456 double* DL_;
457 double* D_;
458 double* DU_;
459 double* DU2_;
460
461 bool UseTranspose_;
462};
463
464// inlined definitions of op() and op[]
465//=========================================================================
466inline double& Ifpack_SerialTriDiMatrix::operator () (int RowIndex, int ColIndex) {
467
468 int diff = ColIndex - RowIndex;
469
470#ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
471 if (ColIndex >= N_ || ColIndex < 0)
472 throw ReportError("Column index = " +toString(ColIndex) +
473 " Out of Range 0 - " + toString(N_-1),-2);
474 if (RowIndex >= N_ || RowIndex < 0)
475 throw ReportError("Row index = " +toString(RowIndex) +
476 " Out of Range 0 - " + toString(N_-1),-2);
477
478 if ( diff > 1 || diff < -1 )
479 throw ReportError("Row index = " +toString(RowIndex) + " differs from Col_Index " + toString(ColIndex) +
480 " Out of Range -1 to 1",-2);
481#endif
482
483 switch (diff) {
484 case -1:
485 // DL
486 return DL_[ColIndex];
487 // break; // unreachable
488 case 0:
489 return D_[ColIndex];
490 // break; // unreachable
491 case 1:
492 return DU_[RowIndex];
493 // break; // unreachable
494 default:
495 throw ReportError("Row index = " +toString(RowIndex) + " differs from Col_Index " + toString(ColIndex) +" Out of Range -1 to 1",1);
496 // return D_[0]; // unreachable
497 }
498 //throw ReportError("Row index = " +toString(RowIndex) + " differs from Col_Index " + toString(ColIndex) + " Out of Range -1 to 1",1); // unreachable
499 // return D_[0]; // unreachable
500}
501//=========================================================================
502inline const double& Ifpack_SerialTriDiMatrix::operator () (int RowIndex, int ColIndex) const {
503 int diff = ColIndex - RowIndex;
504
505#ifdef HAVE_EPETRA_ARRAY_BOUNDS_CHECK
506 if (ColIndex >= N_ || ColIndex < 0)
507 throw ReportError("Column index = " +toString(ColIndex) +
508 " Out of Range 0 - " + toString(N_-1),-2);
509 if (RowIndex >= N_ || RowIndex < 0)
510 throw ReportError("Row index = " +toString(RowIndex) +
511 " Out of Range 0 - " + toString(N_-1),-2);
512 if ( diff > 1 || diff < -1 )
513 throw ReportError("Row index = " +toString(RowIndex) + " differs from Col_Index " + toString(ColIndex) + " Out of Range -1 to 1",-2);
514#endif
515 switch (diff) {
516 case -1:
517 // DL
518 return DL_[ColIndex];
519 // break; // unreachable
520 case 0:
521 return D_[ColIndex];
522 // break; // unreachable
523 case 1:
524 return DU_[RowIndex];
525 // break; // unreachable
526 default:
527 throw ReportError("Row index = " +toString(RowIndex) + " differs from Col_Index " + toString(ColIndex) + " Out of Range -1 to 1",-2);
528 // return D_[0]; // unreachable
529 }
530 // return D_[0]; // unreachable
531}
532
533
534#endif /* EPETRA_SERIALTRIDIMATRIX_H */
Epetra_DataAccess
virtual const char * Label() const
virtual int ReportError(const std::string Message, int ErrorCode) const
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
virtual const char * Label() const
Returns a character string describing the operator.
virtual double OneNorm() const
Computes the 1-Norm of the this matrix (identical to NormOne() method).
virtual double NormOne() const
Computes the 1-Norm of the this matrix.
int Shape(int NumRowCol)
Set dimensions of a Ifpack_SerialTriDiMatrix object; init values to zero.
int N() const
Returns column dimension of system.
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
double * A() const
Returns pointer to the this matrix.
double * DL()
Returns pointer to the this matrix.
bool operator!=(const Ifpack_SerialTriDiMatrix &rhs) const
Inequality operator.
Epetra_DataAccess CV() const
Returns the data access mode of the this matrix.
virtual ~Ifpack_SerialTriDiMatrix()
Ifpack_SerialTriDiMatrix destructor.
virtual int RowColDim() const
Returns the column dimension of operator.
double & operator()(int RowIndex, int ColIndex)
Element access function.
virtual int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
int Random()
Column access function.
int Scale(double ScalarA)
Matrix-Vector multiplication, y = A*x, where 'this' == A.
virtual int ApplyInverse(const Ifpack_SerialTriDiMatrix &X, Ifpack_SerialTriDiMatrix &Y)
Returns the result of a Ifpack_SerialTriDiOperator inverse applied to an Ifpack_SerialTriDiMatrix X i...
virtual void Print(std::ostream &os) const
Print service methods; defines behavior of ostream << operator.
bool operator==(const Ifpack_SerialTriDiMatrix &rhs) const
Comparison operator.
Ifpack_SerialTriDiMatrix & operator=(const Ifpack_SerialTriDiMatrix &Source)
Value copy from one matrix to another.
int Multiply(char TransA, char TransB, double ScalarAB, const Ifpack_SerialTriDiMatrix &A, const Ifpack_SerialTriDiMatrix &B, double ScalarThis)
Matrix-Matrix multiplication, this = ScalarThis*this + ScalarAB*A*B.
virtual double InfNorm() const
Computes the Infinity-Norm of the this matrix (identical to NormInf() method).
virtual double NormInf() const
Computes the Infinity-Norm of the this matrix.
Ifpack_SerialTriDiMatrix & operator+=(const Ifpack_SerialTriDiMatrix &Source)
Add one matrix to another.