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 #include "Epetra_OskiMultiVector.h"
50 #include "Epetra_OskiVector.h"
51 #include "Epetra_CrsMatrix.h"
52 #include "Epetra_OskiPermutation.h"
54 #include "Epetra_Comm.h"
55 extern "C" {
56  #include "oski/oski.h"
57 }
58 
59 class Epetra_OskiVector;
61 class Teuchos_ParameterList;
63 
65 
95  public:
97 
98  Epetra_OskiMatrix(const Epetra_OskiMatrix& Source); //not in use for now
100 
102 
116  Epetra_OskiMatrix(const Epetra_CrsMatrix& Source, const Teuchos::ParameterList& List);
117 
119  virtual ~Epetra_OskiMatrix();
121 
123 
124 
144  int ReplaceMyValues(int MyRow,
145  int NumEntries,
146  double* Values,
147  int* Indices);
148 
150 
167  int SumIntoMyValues(int MyRow,
168  int NumEntries,
169  double* Values,
170  int* Indices);
171 
172 
174 
179  int ExtractDiagonalCopy(Epetra_Vector& Diagonal) const;
180 
182 
192  int ReplaceDiagonalValues(const Epetra_OskiVector& Diagonal);
194 
196 
197 
207  int Multiply(bool TransA,
208  const Epetra_Vector& x,
209  Epetra_Vector& y) const;
210 
212 
223  int Multiply(bool TransA,
224  const Epetra_Vector& x,
225  Epetra_Vector& y,
226  double Alpha,
227  double Beta = 0.0) const;
228 
230 
239  int Multiply(bool TransA,
240  const Epetra_MultiVector& X,
241  Epetra_MultiVector& Y) const;
242 
244 
255  int Multiply(bool TransA,
256  const Epetra_MultiVector& X,
258  double Alpha,
259  double Beta = 0.0) const;
260 
262 
273  int Solve(bool Upper, bool TransA, bool UnitDiagonal, const Epetra_Vector& x, Epetra_Vector &y) const;
274 
276 
286  int Solve(bool TransA, const Epetra_Vector& x, Epetra_Vector& y, double Alpha = 1.0) const;
287 
289 
300  int Solve(bool Upper, bool TransA, bool UnitDiagonal, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
301 
303 
313  int Solve(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y, double Alpha = 1.0) const;
314 
315 
317 
334  int MatTransMatMultiply(bool ATA,
335  const Epetra_Vector& x,
336  Epetra_Vector& y,
337  Epetra_Vector* t,
338  double Alpha = 1.0,
339  double Beta = 0.0) const;
340 
342 
359  int MatTransMatMultiply(bool ATA,
360  const Epetra_MultiVector& X,
363  double Alpha = 1.0,
364  double Beta = 0.0) const;
365 
367 
383  int MultiplyAndMatTransMultiply(bool TransA,
384  const Epetra_Vector& x,
385  Epetra_Vector& y,
386  const Epetra_Vector& w,
387  Epetra_Vector& z,
388  double Alpha = 1.0,
389  double Beta = 0.0,
390  double Omega = 1.0,
391  double Zeta = 0.0) const;
392 
394 
410  int MultiplyAndMatTransMultiply(bool TransA,
411  const Epetra_MultiVector& X,
413  const Epetra_MultiVector& W,
415  double Alpha = 1.0,
416  double Beta = 0.0,
417  double Omega = 1.0,
418  double Zeta = 0.0) const;
419 
421 
439  int MatPowMultiply(bool TransA,
440  const Epetra_Vector& x,
441  Epetra_Vector& y,
443  int Power = 2,
444  double Alpha = 1.0,
445  double Beta = 0.0) const;
446 
448 
463  int MatPowMultiply(bool TransA,
464  const Epetra_Vector& x,
465  Epetra_Vector& y,
466  int Power = 2,
467  double Alpha = 1.0,
468  double Beta = 0.0) const;
470 
472 
473 
507  int SetHint(const Teuchos::ParameterList& List);
508 
509 
511 
545  int SetHintMultiply(bool TransA,
546  double Alpha,
547  const Epetra_OskiMultiVector& InVec,
548  double Beta,
549  const Epetra_OskiMultiVector& OutVec,
550  int NumCalls,
551  const Teuchos::ParameterList& List);
552 
554 
583  int SetHintSolve(bool TransA,
584  double Alpha,
585  const Epetra_OskiMultiVector& Vector,
586  int NumCalls,
587  const Teuchos::ParameterList& List);
588 
590 
630  int SetHintMatTransMatMultiply(bool ATA,
631  double Alpha,
632  const Epetra_OskiMultiVector& InVec,
633  double Beta,
634  const Epetra_OskiMultiVector& OutVec,
635  const Epetra_OskiMultiVector& Intermediate,
636  int NumCalls,
637  const Teuchos::ParameterList& List);
638 
640 
687  int SetHintMultiplyAndMatTransMultiply(bool TransA,
688  double Alpha,
689  const Epetra_OskiMultiVector& InVec,
690  double Beta,
691  const Epetra_OskiMultiVector& OutVec,
692  double Omega,
693  const Epetra_OskiMultiVector& InVec2,
694  double Zeta,
695  const Epetra_OskiMultiVector& OutVec2,
696  int NumCalls,
697  const Teuchos::ParameterList& List);
698 
700 
741  int SetHintPowMultiply(bool TransA,
742  double Alpha,
743  const Epetra_OskiMultiVector& InVec,
744  double Beta,
745  const Epetra_OskiMultiVector& OutVec,
746  const Epetra_OskiMultiVector& Intermediate,
747  int Power,
748  int NumCalls,
749  const Teuchos::ParameterList& List);
750 
752 
757  int TuneMatrix();
759 
761 
762  int IsMatrixTransformed() const;
764 
766  const Epetra_OskiMatrix& ViewTransformedMat() const;
767 
770 
773 
775 
780  char* GetMatrixTransforms() const;
781 
783 
791  int ApplyMatrixTransforms(const char* Transforms);
793 
794  protected:
795 
796  private:
798  oski_matrix_t A_tunable_;
800 };
801 #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...