Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_TriDiContainer.h
Go to the documentation of this file.
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef IFPACK_TRIDICONTAINER_H
44 #define IFPACK_TRIDICONTAINER_H
45 
46 #if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
47 #ifdef __GNUC__
48 #warning "The Ifpack package is deprecated"
49 #endif
50 #endif
51 
52 #include "Ifpack_ConfigDefs.h"
53 #include "Ifpack_Container.h"
56 #include "Epetra_IntSerialDenseVector.h" // Is this needed \cbl
57 #include "Epetra_SerialDenseVector.h"
58 class Epetra_RowMatrix;
59 
61 
125 
126 public:
127 
129 
131  Ifpack_TriDiContainer(const int NumRows_in, const int NumVectors_in = 1) :
132  NumRows_(NumRows_in),
133  NumVectors_(NumVectors_in),
137  ComputeFlops_(0.0),
138  ApplyFlops_(0.0),
139  ApplyInverseFlops_(0.0)
140  {}
141 
144  NumRows_(rhs.NumRows()),
145  NumVectors_(rhs.NumVectors()),
148  IsComputed_(rhs.IsComputed())
149  {
150  Matrix_ = rhs.Matrix();
153  LHS_ = rhs.LHS();
154  RHS_ = rhs.RHS();
155  ID_ = rhs.ID();
156  }
157 
160  {}
162 
164 
167  {
168  if (&rhs == this)
169  return(*this);
170 
171  NumRows_ = rhs.NumRows();
172  NumVectors_ = rhs.NumVectors();
173  IsComputed_ = rhs.IsComputed();
175  Matrix_ = rhs.Matrix();
178  LHS_ = rhs.LHS();
179  RHS_ = rhs.RHS();
180  ID_ = rhs.ID();
181 
182  return(*this);
183  }
184 
186 
188 
190  virtual int NumRows() const;
191 
193  virtual int NumVectors() const
194  {
195  return(NumVectors_);
196  }
197 
199  virtual int SetNumVectors(const int NumVectors_in)
200  {
201  if (NumVectors_ == NumVectors_in)
202  return(0);
203 
204  NumVectors_ = NumVectors_in;
207  // zero out vector elements
208  for (int i = 0 ; i < NumRows_ ; ++i)
209  for (int j = 0 ; j < NumVectors_ ; ++j) {
210  LHS_(i,j) = 0.0;
211  RHS_(i,j) = 0.0;
212  }
213  if (NumRows_!=0)
214  {
216  }
217  return(0);
218  }
219 
221  virtual double& LHS(const int i, const int Vector = 0);
222 
224  virtual double& RHS(const int i, const int Vector = 0);
225 
227 
236  virtual int& ID(const int i);
237 
239  virtual int SetMatrixElement(const int row, const int col,
240  const double value);
241 
243  virtual int SetParameters(Teuchos::ParameterList& /* List */)
244  {
245  return(0);
246  }
247 
249  virtual bool IsInitialized() const
250  {
251  return(IsInitialized_);
252  }
253 
255  virtual bool IsComputed() const
256  {
257  return(IsComputed_);
258  }
259 
261  virtual const char* Label() const
262  {
263  return(Label_.c_str());
264  }
265 
267  virtual int SetKeepNonFactoredMatrix(const bool flag)
268  {
269  KeepNonFactoredMatrix_ = flag;
270  return(0);
271  }
272 
274  virtual bool KeepNonFactoredMatrix() const
275  {
276  return(KeepNonFactoredMatrix_);
277  }
278 
280  virtual const Epetra_SerialDenseMatrix& LHS() const
281  {
282  return(LHS_);
283  }
284 
286  virtual const Epetra_SerialDenseMatrix& RHS() const
287  {
288  return(RHS_);
289  }
290 
292  virtual const Ifpack_SerialTriDiMatrix& Matrix() const
293  {
294  return(Matrix_);
295  }
296 
299  {
300  return(NonFactoredMatrix_);
301  }
302 
304  virtual const Epetra_IntSerialDenseVector& ID() const
305  {
306  return(ID_);
307  }
308 
310 
312  virtual int Initialize();
314 
316  virtual int Compute(const Epetra_RowMatrix& Matrix_in);
317 
319  virtual int Apply();
320 
322  virtual int ApplyInverse();
323 
325 
326  virtual double InitializeFlops() const
327  {
328  return(0.0);
329  }
330 
331  virtual double ComputeFlops() const
332  {
333  return(ComputeFlops_);
334  }
335 
336  virtual double ApplyFlops() const
337  {
338  return(ApplyFlops_);
339  }
340 
341  virtual double ApplyInverseFlops() const
342  {
343  return(ApplyInverseFlops_);
344  }
345 
347  virtual std::ostream& Print(std::ostream& os) const;
348 
349 private:
350 
352  virtual int Extract(const Epetra_RowMatrix& Matrix_in);
353 
355  int NumRows_;
377  std::string Label_;
378 
382  double ApplyFlops_;
385 };
386 
387 #endif
virtual int NumRows() const
Returns the number of rows of the matrix and LHS/RHS.
virtual std::ostream & Print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator&lt;&lt;.
std::string Label_
Label for this object.
Ifpack_TriDiContainer(const Ifpack_TriDiContainer &rhs)
Copy constructor.
double ApplyInverseFlops_
Flops in ApplyInverse().
virtual int NumVectors() const
Returns the number of vectors in LHS/RHS.
Ifpack_TriDiContainer: a class to define containers for dense matrices.
virtual const Epetra_IntSerialDenseVector & ID() const
Returns the integer dense vector of IDs.
virtual int SetMatrixElement(const int row, const int col, const double value)
Set the matrix element (row,col) to value.
bool IsInitialized_
If true, the container has been successfully initialized.
bool KeepNonFactoredMatrix_
If true, keeps a copy of the non-factored matrix.
virtual int Extract(const Epetra_RowMatrix &Matrix_in)
Extract the submatrices identified by the ID set int ID().
virtual const Epetra_SerialDenseMatrix & RHS() const
Returns the dense vector containing the RHS.
virtual double InitializeFlops() const
Returns the flops in Initialize().
Epetra_SerialDenseMatrix RHS_
SerialDense vector representing the RHS.
virtual double & LHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of LHS.
Ifpack_SerialTriDiSolver: A class for solving TriDi linear problems.
virtual const char * Label() const
Returns the label of this container.
virtual const Epetra_SerialDenseMatrix & LHS() const
Returns the dense vector containing the LHS.
int NumVectors_
Number of vectors in the container.
virtual int SetParameters(Teuchos::ParameterList &)
Sets all necessary parameters.
virtual int SetKeepNonFactoredMatrix(const bool flag)
If flag is true, keeps a copy of the non-factored matrix.
Ifpack_SerialTriDiSolver Solver_
TriDi solver (solution will be get using LAPACK).
int NumRows_
Number of rows in the container.
virtual int & ID(const int i)
Returns the ID associated to local row i.
virtual const Ifpack_SerialTriDiMatrix & Matrix() const
Returns the dense matrix or its factors.
virtual double & RHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of RHS.
virtual ~Ifpack_TriDiContainer()
Destructor.
virtual double ComputeFlops() const
Returns the flops in Compute().
Epetra_IntSerialDenseVector ID_
Sets of local rows.
virtual double ApplyFlops() const
Returns the flops in Apply().
virtual int SetNumVectors(const int NumVectors_in)
Sets the number of vectors for LHS/RHS.
#define false
virtual const Ifpack_SerialTriDiMatrix & NonFactoredMatrix() const
Returns the non-factored dense matrix (only if stored).
Ifpack_TriDiContainer & operator=(const Ifpack_TriDiContainer &rhs)
Operator=.
int Reshape(int NumRows, int NumCols)
virtual int Apply()
Apply the matrix to RHS, results are stored in LHS.
double ApplyFlops_
Flops in Apply().
virtual double ApplyInverseFlops() const
Returns the flops in ApplyInverse().
virtual int Compute(const Epetra_RowMatrix &Matrix_in)
Finalizes the linear system matrix and prepares for the application of the inverse.
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
virtual int ApplyInverse()
Apply the inverse of the matrix to RHS, results are stored in LHS.
Ifpack_SerialTriDiMatrix Matrix_
TriDi matrix.
virtual int Initialize()
Initialize the container.
Ifpack_SerialTriDiMatrix NonFactoredMatrix_
TriDi matrix, that contains the non-factored matrix.
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
Sets the pointers for left and right hand side vector(s).
Ifpack_Container: a pure virtual class for creating and solving local linear problems.
bool IsComputed_
If true, the container has been successfully computed.
#define IFPACK_CHK_ERR(ifpack_err)
double ComputeFlops_
Flops in Compute().
virtual bool IsInitialized() const
Returns true is the container has been successfully initialized.
virtual bool KeepNonFactoredMatrix() const
Returns KeepNonFactoredMatrix_.
virtual bool IsComputed() const
Returns true is the container has been successfully computed.
Epetra_SerialDenseMatrix LHS_
SerialDense vector representing the LHS.
Ifpack_TriDiContainer(const int NumRows_in, const int NumVectors_in=1)
Default constructor.