IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Ifpack_SerialTriDiSolver.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_SERIALTRIDISOLVER_H
45 #define IFPACK_SERIALTRIDISOLVER_H
47 //class Epetra_SerialDenseVector;
49 
50 #include "Epetra_Object.h"
51 #include "Epetra_CompObject.h"
52 #include "Epetra_BLAS.h"
53 #include "Teuchos_LAPACK.hpp"
54 
55 
57 
129 //=========================================================================
131  public Epetra_CompObject, public Epetra_BLAS,
132  public Epetra_Object {
133  public:
134 
136 
139 
140 
142  virtual ~Ifpack_SerialTriDiSolver();
144 
146 
147 
150 
152 
157 
158 
160 
161 
163 
165  // void FactorWithEquilibration(bool Flag) {Equilibrate_ = Flag; return;};
166 
168  void SolveWithTranspose(bool Flag) {Transpose_ = Flag; if (Flag) TRANS_ = 'T'; else TRANS_ = 'N'; return;};
169 
171  void SolveToRefinedSolution(bool Flag) {RefineSolution_ = Flag; return;};
172 
174 
177  void EstimateSolutionErrors(bool Flag) ;
179 
181 
182 
184 
187  virtual int Factor(void);
188 
190 
193  virtual int Solve(void);
194 
196 
199  virtual int Invert(void);
200 
202 
205  virtual int ApplyRefinement(void);
206 
208 
211  // int UnequilibrateLHS(void);
212 
214 
220  virtual int ReciprocalConditionEstimate(double & Value);
222 
224 
225 
227  bool Transpose() {return(Transpose_);};
228 
230  bool Factored() {return(Factored_);};
231 
233  bool SolutionErrorsEstimated() {return(SolutionErrorsEstimated_);};
234 
236  bool Inverted() {return(Inverted_);};
237 
239  bool ReciprocalConditionEstimated() {return(ReciprocalConditionEstimated_);};
240 
242  bool Solved() {return(Solved_);};
243 
245  bool SolutionRefined() {return(SolutionRefined_);};
247 
249 
250 
252  Ifpack_SerialTriDiMatrix * Matrix() const {return(Matrix_);};
253 
255  Ifpack_SerialTriDiMatrix * FactoredMatrix() const {return(Factor_);};
256 
258  Epetra_SerialDenseMatrix * LHS() const {return(LHS_);};
259 
261  Epetra_SerialDenseMatrix * RHS() const {return(RHS_);};
262 
264  int N() const {return(N_);};
265 
267  double * A() const {return(A_);};
268 
270  int LDA() const {return(LDA_);};
271 
273  double * B() const {return(B_);};
274 
276  int LDB() const {return(LDB_);};
277 
279  int NRHS() const {return(NRHS_);};
280 
282  double * X() const {return(X_);};
283 
285  int LDX() const {return(LDX_);};
286 
288  double * AF() const {return(AF_);};
289 
291  int LDAF() const {return(LDAF_);};
292 
294  int * IPIV() const {return(IPIV_);};
295 
297  double ANORM() const {return(ANORM_);};
298 
300  double RCOND() const {return(RCOND_);};
301 
303 
305  double ROWCND() const {return(ROWCND_);};
306 
308 
310  double COLCND() const {return(COLCND_);};
311 
313  double AMAX() const {return(AMAX_);};
314 
316  double * FERR() const {return(FERR_);};
317 
319  double * BERR() const {return(BERR_);};
320 
322 
324 
325  virtual void Print(std::ostream& os) const;
328  protected:
329 
330  void AllocateWORK() {if (WORK_==0) {LWORK_ = 4*N_; WORK_ = new double[LWORK_];} return;};
331  void AllocateIWORK() {if (IWORK_==0) IWORK_ = new int[N_]; return;};
332  void InitPointers();
333  void DeleteArrays();
334  void ResetMatrix();
335  void ResetVectors();
336 
337  bool Transpose_;
338  bool Factored_;
339  bool EstimateSolutionErrors_;
340  bool SolutionErrorsEstimated_;
341  bool Solved_;
342  bool Inverted_;
343  bool ReciprocalConditionEstimated_;
344  bool RefineSolution_;
345  bool SolutionRefined_;
346 
347  char TRANS_;
348 
349  int N_;
350  int Min_MN_;
351  int NRHS_;
352  int LDA_;
353  int LDAF_;
354  int LDB_;
355  int LDX_;
356  int INFO_;
357  int LWORK_;
358 
359  int * IPIV_;
360  int * IWORK_;
361 
362  double ANORM_;
363  double RCOND_;
364  double ROWCND_;
365  double COLCND_;
366  double AMAX_;
367 
368  Ifpack_SerialTriDiMatrix * Matrix_;
371  Ifpack_SerialTriDiMatrix * Factor_;
372 
373  double * A_;
374  double * FERR_;
375  double * BERR_;
376  double * AF_;
377  double * WORK_;
378 
379  double * B_;
380  double * X_;
381 
382 
383  private:
384  Teuchos::LAPACK<int,double> lapack;
385  // Ifpack_SerialTriDiSolver copy constructor (put here because we don't want user access)
386 
388  Ifpack_SerialTriDiSolver & operator=(const Ifpack_SerialTriDiSolver& Source);
389 };
390 
391 #endif /* EPETRA_SERIALTRIDISOLVER_H */
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
virtual int Invert(void)
Inverts the this matrix.
Ifpack_SerialTriDiSolver()
Default constructor; matrix should be set using SetMatrix(), LHS and RHS set with SetVectors()...
virtual int Solve(void)
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..
int LDB() const
Returns the leading dimension of the RHS.
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed)...
virtual int ApplyRefinement(void)
Apply Iterative Refinement.
bool Solved()
Returns true if the current set of vectors has been solved.
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
virtual int Factor(void)
Computes the in-place LU factorization of the matrix using the LAPACK routine DGETRF.
void SolveWithTranspose(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor...
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement. ...
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
int N() const
Returns column dimension of system.
int SetMatrix(Ifpack_SerialTriDiMatrix &A)
Sets the pointers for coefficient matrix.
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
Ifpack_SerialTriDiSolver: A class for solving TriDi linear problems.
virtual ~Ifpack_SerialTriDiSolver()
Ifpack_SerialTriDiSolver destructor.
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
void EstimateSolutionErrors(bool Flag)
Causes all solves to estimate the forward and backward solution error.
double * X() const
Returns pointer to current solution.
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
Ifpack_SerialTriDiMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
double * B() const
Returns pointer to current RHS.
double * A() const
Returns pointer to the this matrix.
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place)...
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
virtual int ReciprocalConditionEstimate(double &Value)
Unscales the solution vectors if equilibration was used to solve the system.
int LDA() const
Returns the leading dimension of the this matrix.
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
int LDAF() const
Returns the leading dimension of the factored matrix.
bool Transpose()
Returns true if transpose of this matrix has and will be used.
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
int NRHS() const
Returns the number of current right hand sides and solution vectors.
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
Sets the pointers for left and right hand side vector(s).
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
Ifpack_SerialTriDiMatrix * Matrix() const
Returns pointer to current matrix.
virtual void Print(std::ostream &os) const
Print service methods; defines behavior of ostream &lt;&lt; operator.
int LDX() const
Returns the leading dimension of the solution.
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.