Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_SerialDenseSolver.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 #ifndef EPETRA_SERIALDENSESOLVER_H
45 #define EPETRA_SERIALDENSESOLVER_H
46 
47 #if defined(Epetra_SHOW_DEPRECATED_WARNINGS)
48 #ifdef __GNUC__
49 #warning "The Epetra package is deprecated"
50 #endif
51 #endif
52 
53 
55 #include "Epetra_Object.h"
56 #include "Epetra_CompObject.h"
57 #include "Epetra_BLAS.h"
58 #include "Epetra_LAPACK.h"
59 
60 
62 
139 //=========================================================================
140 class EPETRA_LIB_DLL_EXPORT Epetra_SerialDenseSolver :
141  public Epetra_CompObject, public Epetra_BLAS,
142  public Epetra_LAPACK, public Epetra_Object {
143  public:
144 
146 
149 
150 
152  virtual ~Epetra_SerialDenseSolver();
154 
156 
157 
159  int SetMatrix(Epetra_SerialDenseMatrix & A);
160 
162 
165  int SetVectors(Epetra_SerialDenseMatrix & X, Epetra_SerialDenseMatrix & B);
167 
169 
170 
172 
174  void FactorWithEquilibration(bool Flag) {Equilibrate_ = Flag; return;};
175 
177  void SolveWithTranspose(bool Flag) {Transpose_ = Flag; if (Flag) TRANS_ = 'T'; else TRANS_ = 'N'; return;};
178 
180  void SolveToRefinedSolution(bool Flag) {RefineSolution_ = Flag; return;};
181 
183 
186  void EstimateSolutionErrors(bool Flag) ;
188 
190 
191 
193 
196  virtual int Factor(void);
197 
199 
202  virtual int Solve(void);
203 
205 
208  virtual int Invert(void);
209 
211 
214  virtual int ComputeEquilibrateScaling(void);
215 
217 
220  virtual int EquilibrateMatrix(void);
221 
223 
226  int EquilibrateRHS(void);
227 
228 
230 
233  virtual int ApplyRefinement(void);
234 
236 
239  int UnequilibrateLHS(void);
240 
242 
248  virtual int ReciprocalConditionEstimate(double & Value);
250 
252 
253 
255  bool Transpose() {return(Transpose_);};
256 
258  bool Factored() {return(Factored_);};
259 
261  bool A_Equilibrated() {return(A_Equilibrated_);};
262 
264  bool B_Equilibrated() {return(B_Equilibrated_);};
265 
267  virtual bool ShouldEquilibrate() {ComputeEquilibrateScaling(); return(ShouldEquilibrate_);};
268 
270  bool SolutionErrorsEstimated() {return(SolutionErrorsEstimated_);};
271 
273  bool Inverted() {return(Inverted_);};
274 
276  bool ReciprocalConditionEstimated() {return(ReciprocalConditionEstimated_);};
277 
279  bool Solved() {return(Solved_);};
280 
282  bool SolutionRefined() {return(SolutionRefined_);};
284 
286 
287 
289  Epetra_SerialDenseMatrix * Matrix() const {return(Matrix_);};
290 
292  Epetra_SerialDenseMatrix * FactoredMatrix() const {return(Factor_);};
293 
295  Epetra_SerialDenseMatrix * LHS() const {return(LHS_);};
296 
298  Epetra_SerialDenseMatrix * RHS() const {return(RHS_);};
299 
301  int M() const {return(M_);};
302 
304  int N() const {return(N_);};
305 
307  double * A() const {return(A_);};
308 
310  int LDA() const {return(LDA_);};
311 
313  double * B() const {return(B_);};
314 
316  int LDB() const {return(LDB_);};
317 
319  int NRHS() const {return(NRHS_);};
320 
322  double * X() const {return(X_);};
323 
325  int LDX() const {return(LDX_);};
326 
328  double * AF() const {return(AF_);};
329 
331  int LDAF() const {return(LDAF_);};
332 
334  int * IPIV() const {return(IPIV_);};
335 
337  double ANORM() const {return(ANORM_);};
338 
340  double RCOND() const {return(RCOND_);};
341 
343 
345  double ROWCND() const {return(ROWCND_);};
346 
348 
350  double COLCND() const {return(COLCND_);};
351 
353  double AMAX() const {return(AMAX_);};
354 
356  double * FERR() const {return(FERR_);};
357 
359  double * BERR() const {return(BERR_);};
360 
362  double * R() const {return(R_);};
363 
365  double * C() const {return(C_);};
367 
369 
370  virtual void Print(std::ostream& os) const;
373  protected:
374 
375  void AllocateWORK() {if (WORK_==0) {LWORK_ = 4*N_; WORK_ = new double[LWORK_];} return;};
376  void AllocateIWORK() {if (IWORK_==0) IWORK_ = new int[N_]; return;};
377  void InitPointers();
378  void DeleteArrays();
379  void ResetMatrix();
380  void ResetVectors();
381 
382 
388  bool Factored_;
391  bool Solved_;
392  bool Inverted_;
396 
397  char TRANS_;
398 
399  int M_;
400  int N_;
401  int Min_MN_;
402  int NRHS_;
403  int LDA_;
404  int LDAF_;
405  int LDB_;
406  int LDX_;
407  int INFO_;
408  int LWORK_;
409 
410  int * IPIV_;
411  int * IWORK_;
412 
413  double ANORM_;
414  double RCOND_;
415  double ROWCND_;
416  double COLCND_;
417  double AMAX_;
418 
423 
424  double * A_;
425  double * FERR_;
426  double * BERR_;
427  double * AF_;
428  double * WORK_;
429  double * R_;
430  double * C_;
431 
432  double * B_;
433  double * X_;
434 
435 
436  private:
437  // Epetra_SerialDenseSolver copy constructor (put here because we don't want user access)
438 
441 };
442 
443 #endif /* EPETRA_SERIALDENSESOLVER_H */
bool Transpose()
Returns true if transpose of this matrix has and will be used.
double * B() const
Returns pointer to current RHS.
virtual void Print(std::ostream &os) const
Print object to an output stream Print method.
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
bool B_Equilibrated()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
double * C() const
Returns a pointer to the column scale vector used for equilibration.
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
int M() const
Returns row dimension of system.
void SolveWithTranspose(bool Flag)
If Flag is true, causes all subsequent function calls to work with the transpose of this matrix...
double * R() const
Returns a pointer to the row scaling vector used for equilibration.
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
Epetra_CompObject & operator=(const Epetra_CompObject &src)
Epetra_SerialDenseMatrix: A class for constructing and using real double precision general dense matr...
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
Epetra_SerialDenseMatrix * RHS_
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
Epetra_BLAS: The Epetra BLAS Wrapper Class.
Definition: Epetra_BLAS.h:78
int N() const
Returns column dimension of system.
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
Epetra_Object: The base Epetra class.
Definition: Epetra_Object.h:65
virtual bool ShouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system...
Epetra_SerialDenseMatrix * Matrix_
Epetra_CompObject: Functionality and data that is common to all computational classes.
int LDB() const
Returns the leading dimension of the RHS.
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
double * X() const
Returns pointer to current solution.
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...
Epetra_LAPACK: The Epetra LAPACK Wrapper Class.
Definition: Epetra_LAPACK.h:75
int LDAF() const
Returns the leading dimension of the factored matrix.
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement. ...
int LDX() const
Returns the leading dimension of the solution.
Epetra_SerialDenseMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
bool A_Equilibrated()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
Epetra_SerialDenseMatrix * Factor_
void FactorWithEquilibration(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor...
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
Epetra_SerialDenseMatrix * Matrix() const
Returns pointer to current matrix.
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed)...
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
bool Solved()
Returns true if the current set of vectors has been solved.
int Solve(int, TYPE *, TYPE *, TYPE *)
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place)...
int NRHS() const
Returns the number of current right hand sides and solution vectors.
Epetra_SerialDenseMatrix * LHS_
Epetra_SerialDenseSolver: A class for solving dense linear problems.
double * A() const
Returns pointer to the this matrix.