IFPACK  Development
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 //=========================================================================
130 class EPETRA_LIB_DLL_EXPORT Ifpack_SerialTriDiSolver :
131  public Epetra_CompObject, public Epetra_BLAS,
132  public Epetra_Object {
133  public:
134 
136 
139 
140 
142  virtual ~Ifpack_SerialTriDiSolver();
144 
146 
147 
149  int SetMatrix(Ifpack_SerialTriDiMatrix & A);
150 
152 
155  int SetVectors(Epetra_SerialDenseMatrix & X, Epetra_SerialDenseMatrix & B);
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 */
Epetra_Object
Ifpack_SerialTriDiSolver::LHS
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
Definition: Ifpack_SerialTriDiSolver.h:258
Ifpack_SerialTriDiSolver::FactoredMatrix
Ifpack_SerialTriDiMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
Definition: Ifpack_SerialTriDiSolver.h:255
Ifpack_SerialTriDiSolver::LDAF
int LDAF() const
Returns the leading dimension of the factored matrix.
Definition: Ifpack_SerialTriDiSolver.h:291
Ifpack_SerialTriDiSolver::AMAX
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).
Definition: Ifpack_SerialTriDiSolver.h:313
Ifpack_SerialTriDiSolver::N
int N() const
Returns column dimension of system.
Definition: Ifpack_SerialTriDiSolver.h:264
Ifpack_SerialTriDiSolver::ROWCND
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).
Definition: Ifpack_SerialTriDiSolver.h:305
Ifpack_SerialTriDiSolver::Inverted
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
Definition: Ifpack_SerialTriDiSolver.h:236
Ifpack_SerialTriDiSolver::RCOND
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).
Definition: Ifpack_SerialTriDiSolver.h:300
Ifpack_SerialTriDiMatrix
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
Definition: Ifpack_SerialTriDiMatrix.h:106
Ifpack_SerialTriDiSolver::Solved
bool Solved()
Returns true if the current set of vectors has been solved.
Definition: Ifpack_SerialTriDiSolver.h:242
Ifpack_SerialTriDiSolver::FERR
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
Definition: Ifpack_SerialTriDiSolver.h:316
Ifpack_SerialTriDiSolver::LDX
int LDX() const
Returns the leading dimension of the solution.
Definition: Ifpack_SerialTriDiSolver.h:285
Ifpack_SerialTriDiSolver::LDA
int LDA() const
Returns the leading dimension of the this matrix.
Definition: Ifpack_SerialTriDiSolver.h:270
Ifpack_SerialTriDiSolver::A
double * A() const
Returns pointer to the this matrix.
Definition: Ifpack_SerialTriDiSolver.h:267
Ifpack_SerialTriDiSolver::SolveToRefinedSolution
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement.
Definition: Ifpack_SerialTriDiSolver.h:171
Ifpack_SerialTriDiSolver::NRHS
int NRHS() const
Returns the number of current right hand sides and solution vectors.
Definition: Ifpack_SerialTriDiSolver.h:279
Ifpack_SerialTriDiSolver::LDB
int LDB() const
Returns the leading dimension of the RHS.
Definition: Ifpack_SerialTriDiSolver.h:276
Ifpack_SerialTriDiSolver::RHS
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
Definition: Ifpack_SerialTriDiSolver.h:261
Ifpack_SerialTriDiSolver::Transpose
bool Transpose()
Returns true if transpose of this matrix has and will be used.
Definition: Ifpack_SerialTriDiSolver.h:227
Ifpack_SerialTriDiSolver::X
double * X() const
Returns pointer to current solution.
Definition: Ifpack_SerialTriDiSolver.h:282
Epetra_BLAS
Ifpack_SerialTriDiSolver::SolutionErrorsEstimated
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
Definition: Ifpack_SerialTriDiSolver.h:233
Ifpack_SerialTriDiSolver::Factored
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
Definition: Ifpack_SerialTriDiSolver.h:230
Ifpack_SerialTriDiSolver::SolveWithTranspose
void SolveWithTranspose(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor.
Definition: Ifpack_SerialTriDiSolver.h:168
Ifpack_SerialTriDiSolver::BERR
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
Definition: Ifpack_SerialTriDiSolver.h:319
Ifpack_SerialTriDiSolver::ANORM
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
Definition: Ifpack_SerialTriDiSolver.h:297
Ifpack_SerialTriDiSolver::Matrix
Ifpack_SerialTriDiMatrix * Matrix() const
Returns pointer to current matrix.
Definition: Ifpack_SerialTriDiSolver.h:252
Epetra_CompObject
Ifpack_SerialTriDiSolver::SolutionRefined
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
Definition: Ifpack_SerialTriDiSolver.h:245
Ifpack_SerialTriDiSolver::AF
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place).
Definition: Ifpack_SerialTriDiSolver.h:288
Ifpack_SerialTriDiSolver::COLCND
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
Definition: Ifpack_SerialTriDiSolver.h:310
Epetra_SerialDenseMatrix
Ifpack_SerialTriDiSolver
Ifpack_SerialTriDiSolver: A class for solving TriDi linear problems.
Definition: Ifpack_SerialTriDiSolver.h:130
Epetra_Object::Print
virtual void Print(std::ostream &os) const
Ifpack_SerialTriDiSolver::B
double * B() const
Returns pointer to current RHS.
Definition: Ifpack_SerialTriDiSolver.h:273
Ifpack_SerialTriDiSolver::IPIV
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
Definition: Ifpack_SerialTriDiSolver.h:294
Ifpack_SerialTriDiSolver::ReciprocalConditionEstimated
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
Definition: Ifpack_SerialTriDiSolver.h:239