Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_DenseContainer.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_DENSECONTAINER_H
44 #define IFPACK_DENSECONTAINER_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"
54 #include "Epetra_SerialDenseMatrix.h"
55 #include "Epetra_SerialDenseSolver.h"
56 #include "Epetra_IntSerialDenseVector.h"
57 class Epetra_RowMatrix;
58 
60 
124 
125 public:
126 
128 
130  Ifpack_DenseContainer(const int NumRows_in, const int NumVectors_in = 1) :
131  NumRows_(NumRows_in),
132  NumVectors_(NumVectors_in),
136  ComputeFlops_(0.0),
137  ApplyFlops_(0.0),
138  ApplyInverseFlops_(0.0)
139  {}
140 
143  NumRows_(rhs.NumRows()),
144  NumVectors_(rhs.NumVectors()),
147  IsComputed_(rhs.IsComputed())
148  {
149  Matrix_ = rhs.Matrix();
152  LHS_ = rhs.LHS();
153  RHS_ = rhs.RHS();
154  ID_ = rhs.ID();
155  }
156 
159  {}
161 
163 
166  {
167  if (&rhs == this)
168  return(*this);
169 
170  NumRows_ = rhs.NumRows();
171  NumVectors_ = rhs.NumVectors();
172  IsComputed_ = rhs.IsComputed();
174  Matrix_ = rhs.Matrix();
177  LHS_ = rhs.LHS();
178  RHS_ = rhs.RHS();
179  ID_ = rhs.ID();
180 
181  return(*this);
182  }
183 
185 
187 
189  virtual int NumRows() const;
190 
192  virtual int NumVectors() const
193  {
194  return(NumVectors_);
195  }
196 
198  virtual int SetNumVectors(const int NumVectors_in)
199  {
200  if (NumVectors_ == NumVectors_in)
201  return(0);
202 
203  NumVectors_ = NumVectors_in;
206  // zero out vector elements
207  for (int i = 0 ; i < NumRows_ ; ++i)
208  for (int j = 0 ; j < NumVectors_ ; ++j) {
209  LHS_(i,j) = 0.0;
210  RHS_(i,j) = 0.0;
211  }
212  if (NumRows_!=0)
213  {
215  }
216  return(0);
217  }
218 
220  virtual double& LHS(const int i, const int Vector = 0);
221 
223  virtual double& RHS(const int i, const int Vector = 0);
224 
226 
235  virtual int& ID(const int i);
236 
238  virtual int SetMatrixElement(const int row, const int col,
239  const double value);
240 
242  virtual int SetParameters(Teuchos::ParameterList& /* List */)
243  {
244  return(0);
245  }
246 
248  virtual bool IsInitialized() const
249  {
250  return(IsInitialized_);
251  }
252 
254  virtual bool IsComputed() const
255  {
256  return(IsComputed_);
257  }
258 
260  virtual const char* Label() const
261  {
262  return(Label_.c_str());
263  }
264 
266  virtual int SetKeepNonFactoredMatrix(const bool flag)
267  {
268  KeepNonFactoredMatrix_ = flag;
269  return(0);
270  }
271 
273  virtual bool KeepNonFactoredMatrix() const
274  {
275  return(KeepNonFactoredMatrix_);
276  }
277 
279  virtual const Epetra_SerialDenseMatrix& LHS() const
280  {
281  return(LHS_);
282  }
283 
285  virtual const Epetra_SerialDenseMatrix& RHS() const
286  {
287  return(RHS_);
288  }
289 
291  virtual const Epetra_SerialDenseMatrix& Matrix() const
292  {
293  return(Matrix_);
294  }
295 
298  {
299  return(NonFactoredMatrix_);
300  }
301 
303  virtual const Epetra_IntSerialDenseVector& ID() const
304  {
305  return(ID_);
306  }
307 
309 
311  virtual int Initialize();
313 
315  virtual int Compute(const Epetra_RowMatrix& Matrix_in);
316 
318  virtual int Apply();
319 
321  virtual int ApplyInverse();
322 
324 
325  virtual double InitializeFlops() const
326  {
327  return(0.0);
328  }
329 
330  virtual double ComputeFlops() const
331  {
332  return(ComputeFlops_);
333  }
334 
335  virtual double ApplyFlops() const
336  {
337  return(ApplyFlops_);
338  }
339 
340  virtual double ApplyInverseFlops() const
341  {
342  return(ApplyInverseFlops_);
343  }
344 
346  virtual std::ostream& Print(std::ostream& os) const;
347 
348 private:
349 
351  virtual int Extract(const Epetra_RowMatrix& Matrix_in);
352 
354  int NumRows_;
376  std::string Label_;
377 
381  double ApplyFlops_;
384 };
385 
386 #endif
Ifpack_DenseContainer(const int NumRows_in, const int NumVectors_in=1)
Default constructor.
virtual double ComputeFlops() const
Returns the flops in Compute().
Epetra_SerialDenseMatrix NonFactoredMatrix_
Dense matrix, that contains the non-factored matrix.
virtual const Epetra_IntSerialDenseVector & ID() const
Returns the integer dense vector of IDs.
virtual int Compute(const Epetra_RowMatrix &Matrix_in)
Finalizes the linear system matrix and prepares for the application of the inverse.
virtual const char * Label() const
Returns the label of this container.
virtual std::ostream & Print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator&lt;&lt;.
virtual int SetNumVectors(const int NumVectors_in)
Sets the number of vectors for LHS/RHS.
virtual int SetMatrixElement(const int row, const int col, const double value)
Set the matrix element (row,col) to value.
virtual int SetKeepNonFactoredMatrix(const bool flag)
If flag is true, keeps a copy of the non-factored matrix.
virtual const Epetra_SerialDenseMatrix & Matrix() const
Returns the dense matrix or its factors.
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
virtual int NumVectors() const
Returns the number of vectors in LHS/RHS.
Epetra_SerialDenseMatrix RHS_
Dense vector representing the RHS.
virtual double InitializeFlops() const
Returns the flops in Initialize().
virtual double ApplyInverseFlops() const
Returns the flops in ApplyInverse().
Epetra_SerialDenseSolver Solver_
Dense solver (solution will be get using LAPACK).
virtual bool IsInitialized() const
Returns true is the container has been successfully initialized.
virtual ~Ifpack_DenseContainer()
Destructor.
virtual int Extract(const Epetra_RowMatrix &Matrix_in)
Extract the submatrices identified by the ID set int ID().
bool IsInitialized_
If true, the container has been successfully initialized.
Ifpack_DenseContainer(const Ifpack_DenseContainer &rhs)
Copy constructor.
virtual int SetParameters(Teuchos::ParameterList &)
Sets all necessary parameters.
double ApplyFlops_
Flops in Apply().
virtual double & LHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of LHS.
virtual int & ID(const int i)
Returns the ID associated to local row i.
virtual int Apply()
Apply the matrix to RHS, results are stored in LHS.
virtual bool IsComputed() const
Returns true is the container has been successfully computed.
int NumRows_
Number of rows in the container.
virtual const Epetra_SerialDenseMatrix & LHS() const
Returns the dense vector containing the LHS.
std::string Label_
Label for this object.
Ifpack_DenseContainer & operator=(const Ifpack_DenseContainer &rhs)
Operator=.
Epetra_IntSerialDenseVector ID_
Sets of local rows.
virtual int Initialize()
Initialize the container.
#define false
int Reshape(int NumRows, int NumCols)
virtual int ApplyInverse()
Apply the inverse of the matrix to RHS, results are stored in LHS.
virtual int NumRows() const
Returns the number of rows of the matrix and LHS/RHS.
Epetra_SerialDenseMatrix Matrix_
Dense matrix.
int NumVectors_
Number of vectors in the container.
bool IsComputed_
If true, the container has been successfully computed.
Ifpack_DenseContainer: a class to define containers for dense matrices.
double ComputeFlops_
Flops in Compute().
Ifpack_Container: a pure virtual class for creating and solving local linear problems.
virtual double & RHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of RHS.
double ApplyInverseFlops_
Flops in ApplyInverse().
virtual double ApplyFlops() const
Returns the flops in Apply().
virtual const Epetra_SerialDenseMatrix & NonFactoredMatrix() const
Returns the non-factored dense matrix (only if stored).
#define IFPACK_CHK_ERR(ifpack_err)
bool KeepNonFactoredMatrix_
If true, keeps a copy of the non-factored matrix.
virtual bool KeepNonFactoredMatrix() const
Returns KeepNonFactoredMatrix_.
Epetra_SerialDenseMatrix LHS_
Dense vector representing the LHS.
virtual const Epetra_SerialDenseMatrix & RHS() const
Returns the dense vector containing the RHS.