Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_OskiMatrix.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 // Author: Ian Karlin ikarlin@sandia.gov 05-28-2008
45 
46 #ifndef EPETRA_OSKIMATRIX_H
47 #define EPETRA_OSKIMATRIX_H
48 
49 #if defined(Epetra_SHOW_DEPRECATED_WARNINGS)
50 #ifdef __GNUC__
51 #warning "The Epetra package is deprecated"
52 #endif
53 #endif
54 
55 
56 
57 #include "Epetra_OskiMultiVector.h"
58 #include "Epetra_OskiVector.h"
59 #include "Epetra_CrsMatrix.h"
60 #include "Epetra_OskiPermutation.h"
62 #include "Epetra_Comm.h"
63 extern "C" {
64  #include "oski/oski.h"
65 }
66 
67 class Epetra_OskiVector;
69 class Teuchos_ParameterList;
71 
73 
103  public:
105 
106  Epetra_OskiMatrix(const Epetra_OskiMatrix& Source); //not in use for now
108 
110 
124  Epetra_OskiMatrix(const Epetra_CrsMatrix& Source, const Teuchos::ParameterList& List);
125 
127  virtual ~Epetra_OskiMatrix();
129 
131 
132 
152  int ReplaceMyValues(int MyRow,
153  int NumEntries,
154  double* Values,
155  int* Indices);
156 
158 
175  int SumIntoMyValues(int MyRow,
176  int NumEntries,
177  double* Values,
178  int* Indices);
179 
180 
182 
187  int ExtractDiagonalCopy(Epetra_Vector& Diagonal) const;
188 
190 
200  int ReplaceDiagonalValues(const Epetra_OskiVector& Diagonal);
202 
204 
205 
215  int Multiply(bool TransA,
216  const Epetra_Vector& x,
217  Epetra_Vector& y) const;
218 
220 
231  int Multiply(bool TransA,
232  const Epetra_Vector& x,
233  Epetra_Vector& y,
234  double Alpha,
235  double Beta = 0.0) const;
236 
238 
247  int Multiply(bool TransA,
248  const Epetra_MultiVector& X,
249  Epetra_MultiVector& Y) const;
250 
252 
263  int Multiply(bool TransA,
264  const Epetra_MultiVector& X,
266  double Alpha,
267  double Beta = 0.0) const;
268 
270 
281  int Solve(bool Upper, bool TransA, bool UnitDiagonal, const Epetra_Vector& x, Epetra_Vector &y) const;
282 
284 
294  int Solve(bool TransA, const Epetra_Vector& x, Epetra_Vector& y, double Alpha = 1.0) const;
295 
297 
308  int Solve(bool Upper, bool TransA, bool UnitDiagonal, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
309 
311 
321  int Solve(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y, double Alpha = 1.0) const;
322 
323 
325 
342  int MatTransMatMultiply(bool ATA,
343  const Epetra_Vector& x,
344  Epetra_Vector& y,
345  Epetra_Vector* t,
346  double Alpha = 1.0,
347  double Beta = 0.0) const;
348 
350 
367  int MatTransMatMultiply(bool ATA,
368  const Epetra_MultiVector& X,
371  double Alpha = 1.0,
372  double Beta = 0.0) const;
373 
375 
391  int MultiplyAndMatTransMultiply(bool TransA,
392  const Epetra_Vector& x,
393  Epetra_Vector& y,
394  const Epetra_Vector& w,
395  Epetra_Vector& z,
396  double Alpha = 1.0,
397  double Beta = 0.0,
398  double Omega = 1.0,
399  double Zeta = 0.0) const;
400 
402 
418  int MultiplyAndMatTransMultiply(bool TransA,
419  const Epetra_MultiVector& X,
421  const Epetra_MultiVector& W,
423  double Alpha = 1.0,
424  double Beta = 0.0,
425  double Omega = 1.0,
426  double Zeta = 0.0) const;
427 
429 
447  int MatPowMultiply(bool TransA,
448  const Epetra_Vector& x,
449  Epetra_Vector& y,
451  int Power = 2,
452  double Alpha = 1.0,
453  double Beta = 0.0) const;
454 
456 
471  int MatPowMultiply(bool TransA,
472  const Epetra_Vector& x,
473  Epetra_Vector& y,
474  int Power = 2,
475  double Alpha = 1.0,
476  double Beta = 0.0) const;
478 
480 
481 
515  int SetHint(const Teuchos::ParameterList& List);
516 
517 
519 
553  int SetHintMultiply(bool TransA,
554  double Alpha,
555  const Epetra_OskiMultiVector& InVec,
556  double Beta,
557  const Epetra_OskiMultiVector& OutVec,
558  int NumCalls,
559  const Teuchos::ParameterList& List);
560 
562 
591  int SetHintSolve(bool TransA,
592  double Alpha,
593  const Epetra_OskiMultiVector& Vector,
594  int NumCalls,
595  const Teuchos::ParameterList& List);
596 
598 
638  int SetHintMatTransMatMultiply(bool ATA,
639  double Alpha,
640  const Epetra_OskiMultiVector& InVec,
641  double Beta,
642  const Epetra_OskiMultiVector& OutVec,
643  const Epetra_OskiMultiVector& Intermediate,
644  int NumCalls,
645  const Teuchos::ParameterList& List);
646 
648 
695  int SetHintMultiplyAndMatTransMultiply(bool TransA,
696  double Alpha,
697  const Epetra_OskiMultiVector& InVec,
698  double Beta,
699  const Epetra_OskiMultiVector& OutVec,
700  double Omega,
701  const Epetra_OskiMultiVector& InVec2,
702  double Zeta,
703  const Epetra_OskiMultiVector& OutVec2,
704  int NumCalls,
705  const Teuchos::ParameterList& List);
706 
708 
749  int SetHintPowMultiply(bool TransA,
750  double Alpha,
751  const Epetra_OskiMultiVector& InVec,
752  double Beta,
753  const Epetra_OskiMultiVector& OutVec,
754  const Epetra_OskiMultiVector& Intermediate,
755  int Power,
756  int NumCalls,
757  const Teuchos::ParameterList& List);
758 
760 
765  int TuneMatrix();
767 
769 
770  int IsMatrixTransformed() const;
772 
774  const Epetra_OskiMatrix& ViewTransformedMat() const;
775 
778 
781 
783 
788  char* GetMatrixTransforms() const;
789 
791 
799  int ApplyMatrixTransforms(const char* Transforms);
801 
802  protected:
803 
804  private:
806  oski_matrix_t A_tunable_;
808 };
809 #endif /* EPETRA_OSKIMATRIX_H */
Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in p...
int MatTransMatMultiply(bool ATA, const Epetra_Vector &x, Epetra_Vector &y, Epetra_Vector *t, double Alpha=1.0, double Beta=0.0) const
Performs two matrix vector multiplies of y = Alpha*this^TransA*this*x + Beta*y or y = Alpha*this*this...
int IsMatrixTransformed() const
Returns 1 if the matrix has been reordered by tuning, and 0 if it has not been.
Epetra_OskiPermutation: A class for storing the permutation performed on a Epetra_OskiMatrix.
Epetra_OskiMultiVector: A class for constructing and using dense Oski multi-vectors on a single proce...
const Epetra_OskiPermutation & ViewRowPermutation() const
Returns a read only row/left permutation of the Matrix.
int ExtractDiagonalCopy(Epetra_Vector &Diagonal) const
Returns a copy of the main diagonal in a user-provided vector.
int SetHintMatTransMatMultiply(bool ATA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, const Epetra_OskiMultiVector &Intermediate, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a two matrix-vector multiplies that are composed used by OskiTuneMat to ...
Epetra_OskiVector: A class for constructing and using dense OSKI vectors on a single processor or a s...
oski_matrix_t A_tunable_
char * GetMatrixTransforms() const
Returns a string holding the transformations performed on the matrix when it was tuned.
double ** Values() const
Epetra_Vector: A class for constructing and using dense vectors on a parallel computer.
const Epetra_CrsMatrix * Epetra_View_
int ReplaceDiagonalValues(const Epetra_OskiVector &Diagonal)
Replaces diagonal values of the matrix with those in the user-provided vector.
Epetra_OskiMatrix: A class for constructing and using OSKI Matrices within Epetra. For information on known issues with OSKI see the detailed description.
int SetHintMultiply(bool TransA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a matrix-vector multiply used by OskiTuneMat to optimize the data struct...
int MatPowMultiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y, Epetra_MultiVector &T, int Power=2, double Alpha=1.0, double Beta=0.0) const
Performs a matrix vector multiply of y = Alpha*(this^TransA)^Power*x + Beta*y. This is not implemente...
int TuneMatrix()
Tunes the matrix multiply if its deemed profitable.
int MultiplyAndMatTransMultiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y, const Epetra_Vector &w, Epetra_Vector &z, double Alpha=1.0, double Beta=0.0, double Omega=1.0, double Zeta=0.0) const
Performs the two matrix vector multiplies of y = Alpha*this*x + Beta*y and z = Omega*this^TransA*w + ...
Epetra_OskiMatrix(const Epetra_OskiMatrix &Source)
Copy constructor.
int SetHintMultiplyAndMatTransMultiply(bool TransA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, double Omega, const Epetra_OskiMultiVector &InVec2, double Zeta, const Epetra_OskiMultiVector &OutVec2, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing two matrix-vector multiplies used by OskiTuneMat to optimize the data st...
int Solve(bool Upper, bool TransA, bool UnitDiagonal, const Epetra_Vector &x, Epetra_Vector &y) const
Performs a triangular solve of y = (this^TransA)^-1*x where this is a triangular matrix.
int SetHint(const Teuchos::ParameterList &List)
Stores the hints in List in the matrix structure.
virtual ~Epetra_OskiMatrix()
Destructor.
int SetHintSolve(bool TransA, double Alpha, const Epetra_OskiMultiVector &Vector, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a triangular solve used by OskiTuneMat to optimize the data structure st...
int SumIntoMyValues(int MyRow, int NumEntries, double *Values, int *Indices)
Add this list of entries to existing values for a given local row of the matrix. WARNING: this could ...
const Epetra_OskiPermutation & ViewColumnPermutation() const
Returns a read only column/right permutation of the Matrix.
int SetHintPowMultiply(bool TransA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, const Epetra_OskiMultiVector &Intermediate, int Power, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a matrix-vector multiply performed Power times used by OskiTuneMat to op...
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed r...
int Multiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
Performs a matrix vector multiply of y = this^TransA*x.
int ReplaceMyValues(int MyRow, int NumEntries, double *Values, int *Indices)
Replace current values with this list of entries for a given local row of the matrix. Warning this could be expensive.
int ApplyMatrixTransforms(const char *Transforms)
Replaces the current data structure of the matrix with the one specified in Transforms.
const Epetra_OskiMatrix & ViewTransformedMat() const
Returns the transformed version of InMat if InMat has been transformed. If InMat has not been transfo...